Permutation template

Template for permutations of interval exchange transformations

This file define high level operations on permutations (alphabet, the different rauzy moves, …) shared by reduced and labeled permutations.

AUTHORS:

  • Vincent Delecroix (2008-12-20): initial version

  • Vincent Delecroix (2010-02-11): datatype simplification

TODO:

  • disallow access to stratum, stratum component for permutations with flip

  • construct dynamic Rauzy graphs and paths

  • construct coherent _repr_

surface_dynamics.interval_exchanges.template.FlippedPermutation

alias of Permutation

class surface_dynamics.interval_exchanges.template.FlippedPermutationIET(intervals=None, alphabet=None, reduced=False, flips=None)[source]

Bases: PermutationIET

Template for flipped Abelian permutations.

Warning

Internal class! Do not use directly!

backward_rauzy_move(winner, side=-1)[source]

Returns the permutation before a Rauzy move.

rauzy_move(winner, side='right', inplace=False)[source]

Returns the permutation after a Rauzy move.

class surface_dynamics.interval_exchanges.template.FlippedPermutationLI(intervals=None, alphabet=None, reduced=False, flips=None)[source]

Bases: PermutationLI

Template for flipped quadratic permutations.

Warning

Internal class! Do not use directly!

AUTHORS:

  • Vincent Delecroix (2008-12-20): initial version

backward_rauzy_move(winner, side=-1)[source]

Rauzy move

rauzy_move(winner, side='right', inplace=False)[source]

Rauzy move

class surface_dynamics.interval_exchanges.template.FlippedRauzyDiagram(p, right_induction=True, left_induction=False, left_right_inversion=False, top_bottom_inversion=False, symmetric=False)[source]

Bases: RauzyDiagram

Template for flipped Rauzy diagrams.

AUTHORS:

  • Vincent Delecroix (2009-09-29): initial version

complete(p, reducible=False)[source]

Completion of the Rauzy diagram

Add all successors of p for defined operations in edge_types. Could be used for generating non (strongly) connected Rauzy diagrams. Sometimes, for flipped permutations, the maximal connected graph in all permutations is not strongly connected. Finding such components needs to call most than once the .complete() method.

INPUT:

  • p - a permutation

  • reducible - put or not reducible permutations

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c','c b a',flips='a')
sage: d = p.rauzy_diagram()
sage: d
Rauzy diagram with 3 permutations
sage: p = iet.Permutation('a b c','c b a',flips='b')
sage: d.complete(p)
sage: d
Rauzy diagram with 8 permutations
sage: p = iet.Permutation('a b c','c b a',flips='a')
sage: d.complete(p)
sage: d
Rauzy diagram with 8 permutations
class surface_dynamics.interval_exchanges.template.OrientablePermutationIET(intervals=None, alphabet=None, reduced=False, flips=None)[source]

Bases: PermutationIET

Template for permutation of Interval Exchange Transformation.

Warning

Internal class! Do not use directly!

AUTHOR:

  • Vincent Delecroix (2008-12-20): initial version

arf_invariant()[source]

Returns the Arf invariant of the permutation.

To a permutation pi is associated a quadratic form on the field with 2 elements. The Arf invariant is the total invariant of linear equivalence class of quadratic form of given rank.

Let V be a vector space on the field with two elements FF_2. V there are two equivalence classes of non degenerate quadratic forms. A complete invariant for quadratic forms is the Arf invariant.

For non zero degenerate quadratic forms there are three equivalence classes. If B denotes the bilinear form associated to q then the three classes are as follows

  • the restriction of q to ker(B) is non zero

  • the restriction of q to ker(B) is zero and the spin parity of q on the quotient V/ker(B) is 0

  • the restriction of q to ker(B) is zero and the spin parity of q on the quotient V/ker(B) is 1

The function returns respectively None, 0 or 1 depending on the three alternatives above.

EXAMPLES:

sage: from surface_dynamics import *

Permutations from the odd and even component of H(2,2,2):

sage: a = list(range(10))
sage: b1 = [3,2,4,6,5,7,9,8,1,0]
sage: b0 = [6,5,4,3,2,7,9,8,1,0]
sage: p1 = iet.Permutation(a,b1)
sage: p1.arf_invariant()
1
sage: p0 = iet.Permutation(a,b0)
sage: p0.arf_invariant()
0

Permutations from the odd and even component of H(4,4):

sage: a = list(range(11))
sage: b1 = [3,2,5,4,6,8,7,10,9,1,0]
sage: b0 = [5,4,3,2,6,8,7,10,9,1,0]
sage: p1 = iet.Permutation(a,b1)
sage: p1.arf_invariant()
1
sage: p0 = iet.Permutation(a,b0)
sage: p0.arf_invariant()
0
attached_in_degree()[source]

Returns the degree of the singularity at the right of the interval.

OUTPUT:

  • a positive integer

EXAMPLES:

sage: from surface_dynamics import *

sage: p1 = iet.Permutation('a b c d e f g','d c g f e b a')
sage: p2 = iet.Permutation('a b c d e f g','e d c g f b a')
sage: p1.attached_in_degree()
1
sage: p2.attached_in_degree()
3
attached_out_degree()[source]

Returns the degree of the singularity at the left of the interval.

OUTPUT:

  • a positive integer

EXAMPLES:

sage: from surface_dynamics import *

sage: p1 = iet.Permutation('a b c d e f g','d c g f e b a')
sage: p2 = iet.Permutation('a b c d e f g','e d c g f b a')
sage: p1.attached_out_degree()
3
sage: p2.attached_out_degree()
1
backward_rauzy_move(winner, side='right', inplace=False)[source]

Returns the permutation before a Rauzy move.

INPUT:

  • winner - ‘top’ or ‘bottom’ interval

  • side - ‘right’ or ‘left’ (default: ‘right’) corresponding to the side on which the Rauzy move must be performed.

  • inplace - (default False) whether the Rauzy move is performed inplace (to be used with care since permutations are hashable, set to True if you are sure to know what you are doing)

OUTPUT:

  • a permutation

decompose()[source]

Returns the decomposition as a concatenation of irreducible permutations.

OUTPUT:

a list of permutations

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c','c b a').decompose()[0]
sage: p
a b c
c b a
sage: p1,p2,p3 = iet.Permutation('a b c d e','b a c e d').decompose()
sage: p1
a b
b a
sage: p2
c
c
sage: p3
d e
e d
erase_marked_points()[source]

Returns a permutation equivalent to self but without marked points.

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b','b a')
sage: p.erase_marked_points()
a b
b a
sage: p = iet.Permutation('a b1 b2 c d', 'd c b1 b2 a')
sage: p.erase_marked_points()
a b1 c d
d c b1 a
sage: p = iet.Permutation('a0 a1 b0 b1 c0 c1 d0 d1','d0 d1 c0 c1 b0 b1 a0 a1')
sage: p.erase_marked_points()
a0 b0 c0 d0
d0 c0 b0 a0
sage: p = iet.Permutation('a b y0 y1 x0 x1 c d','c x0 x1 a d y0 y1 b')
sage: p.erase_marked_points()
a b c d
c a d b
sage: p = iet.Permutation('a x y z b','b x y z a')
sage: p.erase_marked_points()
a b
b a
sage: p = iet.Permutation("0 1 2 3 4 5 6","6 0 3 2 4 1 5")
sage: p.stratum()
H_3(4, 0)
sage: p.erase_marked_points().stratum()
H_3(4)
genus()[source]

Returns the genus corresponding to any suspension of self.

The genus can be deduced from the profile (see profile()) p = (p_1,ldots,p_k) of self by the formula: 2g-2 = sum_{i=1}^k (p_i-1).

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c', 'c b a')
sage: p.genus()
1

sage: p = iet.Permutation('a b c d','d c b a')
sage: p.genus()
2
heights_cone(side=None)[source]

Return the cone of heights data.

EXAMPLES:

sage: from surface_dynamics import *
sage: p = iet.Permutation('a b c d', 'd c b a')
sage: C = p.heights_cone()
sage: C
A 4-dimensional polyhedron in QQ^4 defined as the convex hull of 1 vertex and 5 rays
sage: C.rays_list()
[[0, 0, 1, 1], [0, 1, 1, 0], [0, 1, 1, 1], [1, 1, 0, 0], [1, 1, 1, 0]]

sage: p.heights_cone('top').rays_list()
[[0, 0, 1, 1], [0, 1, 1, 0], [1, 1, 0, 0], [1, 1, 1, 0]]
sage: p.heights_cone('bot').rays_list()
[[0, 0, 1, 1], [0, 1, 1, 0], [0, 1, 1, 1], [1, 1, 0, 0]]
intersection_matrix(ring=None)[source]

Returns the intersection matrix.

This d*d antisymmetric matrix is given by the rule :

\[\begin{split}m_{ij} = \begin{cases} 1 & \text{$i < j$ and $\pi(i) > \pi(j)$} \\ -1 & \text{$i > j$ and $\pi(i) < \pi(j)$} \\ 0 & \text{else} \end{cases}\end{split}\]

OUTPUT:

  • a matrix

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c d','d c b a')
sage: p.intersection_matrix()
[ 0  1  1  1]
[-1  0  1  1]
[-1 -1  0  1]
[-1 -1 -1  0]
sage: p = iet.Permutation('1 2 3 4 5','5 3 2 4 1')
sage: p.intersection_matrix()
[ 0  1  1  1  1]
[-1  0  1  0  1]
[-1 -1  0  0  1]
[-1  0  0  0  1]
[-1 -1 -1 -1  0]
sage: p = iet.Permutation('a b c d', 'd c b a')
sage: R = p.rauzy_diagram()
sage: g = R.path(p, *'tbt')
sage: m = g.matrix()
sage: q = g.end()
sage: q.intersection_matrix() == m.transpose() * p.intersection_matrix() * m
True
invariant_density_rauzy(winner=None, var='x')[source]

Return the invariant density for the Rauzy induction.

Goes via the zippered rectangle construction of [Vee1982].

EXAMPLES:

sage: from surface_dynamics import iet
sage: f = iet.Permutation('a b c d', 'd c b a').invariant_density_rauzy()
sage: f
(1)/((x2 + x3)*(x1 + x2)*(x1 + x2 + x3)*(x0 + x1)) + (1)/((x2 + x3)*(x1 + x2)*(x0 + x1)*(x0 + x1 + x2))

sage: f_top = iet.Permutation('a b c d', 'd c b a').invariant_density_rauzy('top')
sage: f_top
(1)/((x2 + x3)*(x1 + x2)*(x0 + x1)*(x0 + x1 + x2))
sage: f_bot = iet.Permutation('a b c d', 'd c b a').invariant_density_rauzy('bot')
sage: f_bot
(1)/((x2 + x3)*(x1 + x2)*(x1 + x2 + x3)*(x0 + x1))

sage: f == f_bot + f_top
True
is_cylindric()[source]

Returns True if the permutation is cylindric

A permutation pi is cylindric if pi(1) = n or pi(n) = 1. The name cylindric comes from geometry. A cylindric permutation has a suspension which is a flat surface with a completely periodic horizontal direction which is made of only one cylinder.

EXAMPLES:

sage: from surface_dynamics import *

sage: iet.Permutation('1 2 3','3 2 1').is_cylindric()
True
sage: iet.Permutation('1 2 3','3 1 2').is_cylindric()
True
sage: iet.Permutation('1 2 3 4','3 1 2 4').is_cylindric()
False
is_hyperelliptic()[source]

Returns True if the permutation is in the class of the symmetric permutations (with eventual marked points).

This is equivalent to say that the suspension lives in an hyperelliptic stratum of Abelian differentials H_hyp(2g-2) or H_hyp(g-1, g-1) with some marked points.

EXAMPLES:

sage: from surface_dynamics import *

sage: iet.Permutation('a b c d','d c b a').is_hyperelliptic()
True
sage: iet.Permutation('0 1 2 3 4 5','5 2 1 4 3 0').is_hyperelliptic()
False
is_identity()[source]

Returns True if self is the identity.

EXAMPLES:

sage: from surface_dynamics import *

sage: iet.Permutation("a b","a b",reduced=False).is_identity()
True
sage: iet.Permutation("a b","a b",reduced=True).is_identity()
True
sage: iet.Permutation("a b","b a",reduced=False).is_identity()
False
sage: iet.Permutation("a b","b a",reduced=True).is_identity()
False
is_standard()[source]

Test if the permutation is standard

A permutation pi is standard if ‘pi(n) = 1` and pi(1) = n.

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c d','d c b a')
sage: p.is_standard()
True
sage: p = p.rauzy_move('top')
sage: p.is_standard()
False
marked_profile()[source]

Returns the marked profile of the permutation

The marked profile of a permutation corresponds to the integer partition associated to the angles of conical singularities in the suspension together with a data associated to the endpoint called marking.

If the left endpoint and the right endpoint of the interval associated to the permutation, then the marking is of type one and consists in a couple (m,a) such that m is the angle of the conical singularity and a is the angle between the outgoing separatrix associated to the left endpoint and the incoming separatrix associated to the right endpoint. A marking of type one is denoted (m|a).

If the left endpoint and the right endpoint are two different conical singularities in the suspension the marking is of type two and consists in a couple (m_l,m_r) where m_l (resp. m_r) is the conical angle of the singularity at the left endpoint (resp. right endpoint). A marking of type two is denoted m_l o m_r

EXAMPLES:

sage: from surface_dynamics import *

The irreducible permutation on 1 interval has marked profile of type 2 with data (0,0):

sage: p = iet.Permutation('a','a')
sage: p.marked_profile()
0o0 []

Permutations in H(3,1) with all possible profiles:

sage: p = iet.Permutation('a b c d e f g','b g a c f e d')
sage: p.interval_diagram()
[['a', ('b', 'a'), ('g', 'd'), 'e', 'f', 'g', 'b', 'c'], ['c', 'd', 'e', 'f']]
sage: p.marked_profile()
4|0 [4, 2]

sage: p = iet.Permutation('a b c d e f g','c a g d f b e')
sage: p.interval_diagram()
[['a', 'b', 'f', 'g'], ['c', 'd', ('g', 'e'), 'f', 'd', 'e', 'b', ('c', 'a')]]
sage: p.marked_profile()
4|1 [4, 2]

sage: p = iet.Permutation('a b c d e f g','e b d g c a f')
sage: p.interval_diagram()
[['a', 'b', 'e', 'f'], ['c', 'd', 'b', 'c', ('g', 'f'), 'g', 'd', ('e', 'a')]]
sage: p.marked_profile()
4|2 [4, 2]

sage: p = iet.Permutation('a b c d e f g', 'e c g b a f d')
sage: p.interval_diagram()
[['a', 'b', ('g', 'd'), ('e', 'a'), 'b', 'c', 'e', 'f'], ['c', 'd', 'f', 'g']]
sage: p.marked_profile()
4|3 [4, 2]

sage: p = iet.Permutation('a b c d e f g', 'f d c a g e b')
sage: p.interval_diagram()
[['a', 'b', 'e', ('f', 'a'), 'c', 'd', 'f', 'g'], ['c', 'd', 'e', ('g', 'b')]]
sage: p.marked_profile()
4o2 [4, 2]
marking()[source]

Return the marking induced by the two sides of the interval

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c d e f','f a e b d c')
sage: p.marking()
5|0
sage: p = iet.Permutation('0 1 2 3 4 5 6','3 2 4 6 5 1 0')
sage: p.marking()
3o3
masur_polygon(lengths, heights)[source]

Return the Masur polygon for the given lengths and heights

EXAMPLES:

sage: from surface_dynamics import iet

sage: p = iet.Permutation('a b c', 'c b a')
sage: S = p.masur_polygon([1,4,2], [2,0,-1])  # optional: sage_flatsurf
sage: stratum = S.stratum()                   # optional: sage_flatsurf # random
sage: stratum                                 # optional: sage_flatsurf
H_1(0^2)
sage: S                                       # optional: sage_flatsurf
Translation Surface in H_1(0^2) built from 2 isosceles triangles and 2 triangles

Generic construction using suspension cone:

sage: p = iet.Permutation('a b c d e f g h i', 'b g f c a d i e h')
sage: x = polygen(QQ)
sage: poly = x^3 - x - 1
sage: emb = AA.polynomial_root(poly, RIF(1.3,1.4))
sage: K = NumberField(poly, 'a', embedding=emb)
sage: a = K.gen()
sage: R = [r.vector() for r in p.suspension_cone().rays()]
sage: C = [1, a, a+1, 2-a, 2, 1, a, a, 1, a-1, 1]
sage: H = sum(c*r for c,r in zip(C,R))
sage: H
(a + 2, -2, 2, -2*a + 2, 3*a, -a - 4, 0, -a + 1, -2*a - 1)
sage: L = [1+a**2, 2*a**2-1, 1, 1, 1+a, a**2, a-1, a-1, 2]
sage: S = p.masur_polygon(L, H)   # optional: sage_flatsurf
sage: S                           # optional: sage_flatsurf
Translation Surface in H_3(1^4) built from 15 triangles and a right triangle
order_of_rauzy_action(winner, side=None)[source]

Returns the order of the action of a Rauzy move.

INPUT:

  • winner - string 'top' or 'bottom'

  • side - string 'left' or 'right'

OUTPUT:

An integer corresponding to the order of the Rauzy action.

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c d','d a c b')
sage: p.order_of_rauzy_action('top', 'right')
3
sage: p.order_of_rauzy_action('bottom', 'right')
2
sage: p.order_of_rauzy_action('top', 'left')
1
sage: p.order_of_rauzy_action('bottom', 'left')
3
profile()[source]

Returns the profile of the permutation

EXAMPLES:

sage: from surface_dynamics import *

sage: iet.Permutation('a b c d','d c b a').profile()
[3]
sage: iet.Permutation('a b c d e','e d c b a').profile()
[2, 2]
rauzy_move(winner, side='right', inplace=False)[source]

Returns the permutation after a Rauzy move.

INPUT:

  • winner - ‘top’ or ‘bottom’ interval

  • side - ‘right’ or ‘left’ (default: ‘right’) corresponding to the side on which the Rauzy move must be performed.

  • inplace - (default False) whether the Rauzy move is performed inplace (to be used with care since permutations are hashable, set to True if you are sure to know what you are doing)

OUTPUT:

  • a permutation

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b','b a')
sage: p.rauzy_move(winner='top', side='right') == p
True
sage: p.rauzy_move(winner='bottom', side='right') == p
True
sage: p.rauzy_move(winner='top', side='left') == p
True
sage: p.rauzy_move(winner='bottom', side='left') == p
True

The options winner can be shortened to ‘t’, ‘b’ and ‘r’, ‘l’. As you can see in the following example:

sage: p = iet.Permutation('a b c','c b a')
sage: p.rauzy_move(winner='t', side='r')
a b c
c a b
sage: p.rauzy_move(winner='b', side='r')
a c b
c b a
sage: p.rauzy_move(winner='t', side='l')
a b c
b c a
sage: p.rauzy_move(winner='b', side='l')
b a c
c b a

This works as well for reduced permutations:

sage: p = iet.Permutation('a b c d','d b c a',reduced=True)
sage: p.rauzy_move('t')
a b c d
d a b c

If Rauzy induction is not well defined, an error is raised:

sage: p = iet.Permutation('a b', 'a b')
sage: p.rauzy_move('t')
Traceback (most recent call last):
...
ValueError: Rauzy induction is not well defined

Test the inplace option:

sage: p = iet.Permutation('a b c d', 'd c b a')
sage: q = p.rauzy_move('t', inplace=True)
sage: assert q is p
sage: p
a b c d
d a c b
sage: q = p.rauzy_move('b', inplace=True)
sage: assert q is p
stratum()[source]

Returns the strata in which any suspension of this permutation lives.

OUTPUT:

  • a stratum of Abelian differentials

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c', 'c b a')
sage: p.stratum()
H_1(0^2)

sage: p = iet.Permutation('a b c d', 'd a b c')
sage: p.stratum()
H_1(0^3)

sage: p = iet.Permutation(list(range(9)), [8,5,2,7,4,1,6,3,0])
sage: p.stratum()
H_3(1^4)

sage: a = 'a b c d e f g h i j'
sage: b3 = 'd c g f e j i h b a'
sage: b2 = 'd c e g f j i h b a'
sage: b1 = 'e d c g f h j i b a'
sage: p3 = iet.Permutation(a, b3)
sage: p3.stratum()
H_4(3, 2, 1)
sage: p2 = iet.Permutation(a, b2)
sage: p2.stratum()
H_4(3, 2, 1)
sage: p1 = iet.Permutation(a, b1)
sage: p1.stratum()
H_4(3, 2, 1)

AUTHORS:

  • Vincent Delecroix (2008-12-20)

stratum_component()[source]

Returns a connected components of a stratum.

EXAMPLES:

sage: from surface_dynamics import *

Permutations from the stratum H(6):

sage: a = list(range(8))
sage: b_hyp = [7,6,5,4,3,2,1,0]
sage: b_odd = [3,2,5,4,7,6,1,0]
sage: b_even = [5,4,3,2,7,6,1,0]
sage: p_hyp = iet.Permutation(a, b_hyp)
sage: p_odd = iet.Permutation(a, b_odd)
sage: p_even = iet.Permutation(a, b_even)
sage: p_hyp.stratum_component()
H_4(6)^hyp
sage: p_odd.stratum_component()
H_4(6)^odd
sage: p_even.stratum_component()
H_4(6)^even

Permutations from the stratum H(4,4):

sage: a = list(range(11))
sage: b_hyp = [10,9,8,7,6,5,4,3,2,1,0]
sage: b_odd = [3,2,5,4,6,8,7,10,9,1,0]
sage: b_even = [5,4,3,2,6,8,7,10,9,1,0]
sage: p_hyp = iet.Permutation(a,b_hyp)
sage: p_odd = iet.Permutation(a,b_odd)
sage: p_even = iet.Permutation(a,b_even)
sage: p_hyp.stratum() == Stratum([4,4], k=1)
True
sage: p_hyp.stratum_component()
H_5(4^2)^hyp
sage: p_odd.stratum() == Stratum([4,4], k=1)
True
sage: p_odd.stratum_component()
H_5(4^2)^odd
sage: p_even.stratum() == Stratum([4,4], k=1)
True
sage: p_even.stratum_component()
H_5(4^2)^even

As for stratum you can specify that you want to attach the singularity on the left of the interval using the option marked_separatrix:

sage: a = list(range(1,10))
sage: b_odd = [4,3,6,5,7,9,8,2,1]
sage: b_even = [6,5,4,3,7,9,8,2,1]
sage: p_odd = iet.Permutation(a,b_odd)
sage: p_even = iet.Permutation(a,b_even)
sage: p_odd.stratum_component()
H_4(4, 2)^odd
sage: p_even.stratum_component()
H_4(4, 2)^even
suspension_cone(winner=None)[source]

Return the cone of suspension data.

A suspension data tau for a permutation (pi_{top}, pi_{bot}) on the alphabet mathcal{A} is a real vector in RR^mathcal{A} so that

\[\forall 1 \leq k < d,\, \sum_{\beta: \pi_{top}(\beta) \leq k} \tau_\beta > 0 \quad \text{and} \quad \sum_{\beta: \pi_{bot}(\beta) \leq k} \tau_\beta < 0.\]

A suspension data determines half of a zippered rectangle construction. The other half is the length data that is a positive vector in RR^mathcal{A}.

INPUT:

  • winner - (optional) either None, "top" or "bottom". If not None , then return only half of the suspension cone corresponding to data that either comes from a top or bottom Rauzy induction.

See also

heights_cone()

EXAMPLES:

sage: from surface_dynamics import *
sage: p = iet.Permutation('a b c d e f', 'e c b f d a')
sage: H = p.suspension_cone()
sage: H.dimension()
6
sage: rays = [r.vector() for r in H.rays()]
sage: r = sum(randint(1,5)*ray for ray in rays)
sage: r[0]>0 and r[0]+r[1] > 0 and r[0]+r[1]+r[2] > 0
True
sage: r[0]+r[1]+r[2]+r[3]>0
True
sage: r[0]+r[1]+r[2]+r[3]+r[4]>0
True
sage: r[4]<0 and r[4]+r[2]<0 and r[4]+r[2]+r[1] < 0
True
sage: r[4]+r[2]+r[1]+r[5]<0
True
sage: r[4]+r[2]+r[1]+r[5]+r[3]<0
True

The construction also works with reduced permutations (ie not carrying labels):

sage: p = iet.Permutation('a b c d', 'd c b a', reduced=True)
sage: H = p.suspension_cone()
sage: r = sum(r.vector() for r in H.rays())
sage: r[0] > 0 and r[0]+r[1] > 0 and r[0]+r[1]+r[2] > 0
True
sage: r[3] < 0 and r[3]+r[2] < 0 and r[3]+r[2]+r[1] < 0
True
to_cylindric()[source]

Returns a cylindric permutation in the same Rauzy class.

A permutation is cylindric if the first letter in the top interval is also the last letter of the bottom interval or if the last letter of the top interval is the first letter of the bottom interval.

to_origami()[source]

Return the origami associated to a cylindric permutation.

EXAMPLES:

sage: from surface_dynamics import iet
sage: p = iet.Permutation('a b', 'b a')
sage: p.to_origami()
(1)
(1)

sage: p = iet.Permutation('a b c e d f g', 'f e b g d c a')
sage: p.to_origami()
(1,2,3,4,5,6)
(1,3,2,6,4,5)
sage: assert p.stratum_component() == p.to_origami().stratum_component()
to_permutation()[source]

Returns the permutation as an element of the symmetric group.

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c','c b a')
sage: p.to_permutation()
[3, 2, 1]
sage: p = Permutation([2,4,1,3])
sage: q = iet.Permutation(p)
sage: q.to_permutation() == p
True
to_standard()[source]

Returns a standard permutation in the same Rauzy class.

class surface_dynamics.interval_exchanges.template.OrientablePermutationLI(intervals=None, alphabet=None, reduced=False, flips=None)[source]

Bases: PermutationLI

Template for quadratic permutation.

Warning

Internal class! Do not use directly!

AUTHOR:

  • Vincent Delecroix (2008-12-20): initial version

backward_rauzy_move(winner, side='top')[source]

Return the permutation before the Rauzy move.

rauzy_move(winner, side='right', inplace=False)[source]

Returns the permutation after a Rauzy move.

stratum()[source]

Returns the stratum associated to self

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.GeneralizedPermutation('a b b','c c a')
sage: p.stratum()
Q_0(-1^4)
class surface_dynamics.interval_exchanges.template.Permutation(intervals=None, alphabet=None, reduced=False, flips=None)[source]

Bases: SageObject

Template for all permutations.

Warning

Internal class! Do not use directly!

This class implement generic algorithm (stratum, connected component, …) and unfies all its children.

It has four attributes

  • _alphabet – the alphabet on which the permutation is defined. Be careful, it might have a different cardinality as the size of the permutation!

  • _twin – the permutation

  • _labels – None or the list of labels

  • _flips – None or the list of flips (each flip is either 1 or -1)

The datatype for _twin differs for IET and LI (TODO: unify).

alphabet(data=None)[source]

Manages the alphabet of self.

If there is no argument, the method returns the alphabet used. If the argument could be converted to an alphabet, this alphabet will be used.

INPUT:

  • data - None or something that could be converted to an alphabet

OUTPUT:

  • either None or the current alphabet

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b','a b')
sage: p.alphabet([0,1])
sage: p.alphabet() == Alphabet([0,1])
True
sage: p
0 1
0 1
sage: p.alphabet("cd")
sage: p.alphabet() == Alphabet(['c','d'])
True
sage: p
c d
c d
cover(perms, as_tuple=False)[source]

Return a covering of this permutation.

INPUT:

  • perms - a list of permutations that describe the gluings

  • as_tuple - whether permutations need to be considered as 1-based (default) or 0-based.

EXAMPLES:

sage: from surface_dynamics import iet
sage: p = iet.Permutation('a b', 'b a')
sage: p.cover(['(1,2)', '(1,3)'])
Covering of degree 3 of the permutation:
a b
b a

sage: p.cover([[1,0,2], [2,1,0]], as_tuple=True)
Covering of degree 3 of the permutation:
a b
b a

sage: p = iet.GeneralizedPermutation('a a b b','c c')
sage: q = p.cover(['(0,1)', [], [2,1,0]], as_tuple=True)
sage: q
Covering of degree 3 of the permutation:
a a b b
c c
sage: q.covering_data('a')
(1,2)
sage: q.covering_data('b')
()
sage: q.covering_data('c')
(1,3)
flips()[source]

Returns the list of flips.

If the permutation is not a flipped permutations then None is returned.

EXAMPLES:

sage: from surface_dynamics import *

sage: iet.Permutation('a b c', 'c b a').flips()
[]
sage: iet.Permutation('a b c', 'c b a', flips='ac').flips()
['a', 'c']
sage: iet.GeneralizedPermutation('a a', 'b b', flips='a').flips()
['a']
sage: iet.GeneralizedPermutation('a a','b b', flips='b', reduced=True).flips()
['b']
horizontal_inverse(inplace=False)

Return the top-bottom inverse.

You can use also use the shorter .tb_inverse().

There are two other symmetries of permutation which are accessible via the methods Permutation.left_right_inverse() and Permutation.symmetric().

OUTPUT: a permutation

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b','b a')
sage: p.top_bottom_inverse()
b a
a b
sage: p = iet.Permutation('a b','b a',reduced=True)
sage: p.top_bottom_inverse() == p
True
sage: p = iet.Permutation('a b c d','c d a b')
sage: p.top_bottom_inverse()
c d a b
a b c d
interval_diagram(glue_ends=True, sign=False)[source]

Return the interval diagram of self.

The result is in random order.

INPUT:

  • glue_ends - bool (default: True)

  • sign - bool (default: False)

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c','c b a')
sage: p.interval_diagram()
[['b', ('c', 'a')], ['b', ('c', 'a')]]

sage: p = iet.Permutation('a b c','c a b')
sage: p.interval_diagram()
[['a', 'b'], [('c', 'b'), ('c', 'a')]]

sage: p = iet.GeneralizedPermutation('a a','b b c c')
sage: p.interval_diagram()
[['a'], [('b', 'a', 'c')], ['b'], ['c']]

sage: p = iet.GeneralizedPermutation('a a b b','c c')
sage: p.interval_diagram()
[['a'], [('b', 'c', 'a')], ['b'], ['c']]
sage: p.interval_diagram(sign=True)
[[('a', -1)], [(('b', 1), ('c', 1), ('a', 1))], [('b', -1)], [('c', -1)]]

sage: p = iet.GeneralizedPermutation((0,1,0,2),(3,2,4,1,4,3))
sage: p.interval_diagram()
[[0, 1, 4, 1], [2, 3, 4, (2, 3, 0)]]
sage: p.interval_diagram(sign=True)
[[(0, -1), (1, 1), (4, -1), (1, -1)],
 [(2, 1), (3, -1), (4, 1), ((2, -1), (3, 1), (0, 1))]]

sage: p = iet.GeneralizedPermutation('a b c d b', 'e d f e a f c', flips='bdf')
sage: p.interval_diagram()
[['a', 'b', 'd', 'e', 'f', 'c', ('b', 'c'), 'd', 'f'], [('e', 'a')]]
left_right_inverse(inplace=False)[source]

Return the left-right inverse.

The left-right inverse of a permutation, is the permutation obtained by reversing the order of the underlying ordering.

You can also use the shorter .lr_inverse()

There are two other symmetries of permutation which are accessible via the methods Permutation.top_bottom_inverse() and Permutation.symmetric().

OUTPUT: a permutation

EXAMPLES:

sage: from surface_dynamics import *

For labelled permutations:

sage: p = iet.Permutation('a b c','c a b')
sage: p.left_right_inverse()
c b a
b a c
sage: p = iet.Permutation('a b c d','c d a b')
sage: p.left_right_inverse()
d c b a
b a d c

for reduced permutations:

sage: p = iet.Permutation('a b c','c a b',reduced=True)
sage: p.left_right_inverse()
a b c
b c a
sage: p = iet.Permutation('a b c d','c d a b',reduced=True)
sage: p.left_right_inverse()
a b c d
c d a b

for labelled quadratic permutations:

sage: p = iet.GeneralizedPermutation('a a','b b c c')
sage: p.left_right_inverse()
a a
c c b b

for reduced quadratic permutations:

sage: p = iet.GeneralizedPermutation('a a','b b c c',reduced=True)
sage: p.left_right_inverse() == p
True
length(interval=None)[source]

Returns the 2-uple of lengths.

p.length() is identical to (p.length_top(), p.length_bottom()) If an interval is specified, it returns the length of the specified interval.

INPUT:

  • interval - None, ‘top’ (or ‘t’ or 0) or ‘bottom’ (or ‘b’ or 1)

OUTPUT:

integer or 2-uple of integers – the corresponding lengths

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c','c b a',reduced=False)
sage: p.length()
(3, 3)
sage: p = iet.Permutation('a b c','c b a',reduced=True)
sage: p.length()
(3, 3)

sage: p = iet.GeneralizedPermutation('a a b','c d c b d',reduced=False)
sage: p.length()
(3, 5)
sage: p = iet.GeneralizedPermutation('a a b','c d c b d',reduced=True)
sage: p.length()
(3, 5)
length_bottom()[source]

Returns the number of intervals in the bottom segment.

OUTPUT:

integer – the length of the bottom segment

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c','c b a',reduced=True)
sage: p.length_bottom()
3
sage: p = iet.Permutation('a b c','c b a',reduced=False)
sage: p.length_bottom()
3
sage: p = iet.GeneralizedPermutation('a a b','c d c b d',reduced=True)
sage: p.length_bottom()
5
sage: p = iet.GeneralizedPermutation('a a b','c d c b d',reduced=False)
sage: p.length_bottom()
5
length_top()[source]

Returns the number of intervals in the top segment.

OUTPUT:

integer – the length of the top segment

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c','c b a', reduced=True)
sage: p.length_top()
3
sage: p = iet.Permutation('a b c','c b a', reduced=False)
sage: p.length_top()
3
sage: p = iet.GeneralizedPermutation('a a b','c d c b d',reduced=True)
sage: p.length_top()
3
sage: p = iet.GeneralizedPermutation('a a b','c d c d b',reduced=False)
sage: p.length_top()
3
sage: p = iet.GeneralizedPermutation('a b c b d c d', 'e a e',reduced=True)
sage: p.length_top()
7
sage: p = iet.GeneralizedPermutation('a b c d b c d', 'e a e', reduced=False)
sage: p.length_top()
7
letters()[source]

Returns the list of letters of the alphabet used for representation.

The letters used are not necessarily the whole alphabet (for example if the alphabet is infinite).

OUTPUT: a list of labels

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation([1,2],[2,1])
sage: p.alphabet(Alphabet(name="NN"))
sage: p
0 1
1 0
sage: p.letters()
[0, 1]

sage: p = iet.GeneralizedPermutation('a a b','b c c')
sage: p.letters()
['a', 'b', 'c']
sage: p.alphabet(range(10))
sage: p.letters()
[0, 1, 2]
sage: p._remove_interval(0, 2)
sage: p.letters()
[0, 2]

sage: p = iet.GeneralizedPermutation('a a b', 'b c c', reduced=True)
sage: p.letters()
['a', 'b', 'c']

For permutations with flips, the letters appear as pairs made of an element of the alphabet and the flip:

sage: p = iet.Permutation('A B C D', 'D C A B', flips='AC')
sage: p.letters()
['A', 'B', 'C',  'D']

sage: p = iet.GeneralizedPermutation('A A B', 'B C C', flips='B')
sage: p.letters()
['A', 'B', 'C']
lr_inverse(inplace=False)

Return the left-right inverse.

The left-right inverse of a permutation, is the permutation obtained by reversing the order of the underlying ordering.

You can also use the shorter .lr_inverse()

There are two other symmetries of permutation which are accessible via the methods Permutation.top_bottom_inverse() and Permutation.symmetric().

OUTPUT: a permutation

EXAMPLES:

sage: from surface_dynamics import *

For labelled permutations:

sage: p = iet.Permutation('a b c','c a b')
sage: p.left_right_inverse()
c b a
b a c
sage: p = iet.Permutation('a b c d','c d a b')
sage: p.left_right_inverse()
d c b a
b a d c

for reduced permutations:

sage: p = iet.Permutation('a b c','c a b',reduced=True)
sage: p.left_right_inverse()
a b c
b c a
sage: p = iet.Permutation('a b c d','c d a b',reduced=True)
sage: p.left_right_inverse()
a b c d
c d a b

for labelled quadratic permutations:

sage: p = iet.GeneralizedPermutation('a a','b b c c')
sage: p.left_right_inverse()
a a
c c b b

for reduced quadratic permutations:

sage: p = iet.GeneralizedPermutation('a a','b b c c',reduced=True)
sage: p.left_right_inverse() == p
True
regular_cover(G, elts)[source]

Return a regular (or normal) cover of this permutation.

EXAMPLES:

sage: from surface_dynamics import iet

sage: p = iet.Permutation('a b c d', 'd c b a')
sage: G = SymmetricGroup(4)
sage: ga = G('(1,2)')
sage: gb = G('(1,3,4)')
sage: gc = G('(1,3)')
sage: gd = G('()')
sage: pp = p.regular_cover(G, [ga, gb, gc, gd])
sage: pp
Regular cover of degree 24 with group Symmetric group of order 4! as a permutation group of the permutation:
a b c d
d c b a
sage: pp.genus()
33

Additive groups are converted to multiplicative groups:

sage: p = iet.GeneralizedPermutation('a a b', 'b c c')
sage: G = Zmod(5)
sage: p.regular_cover(G, [0, 1, 2])
Regular cover of degree 5 with group Multiplicative Abelian group isomorphic to C5 of the permutation:
a a b
b c c
str(sep='\n', align=None)[source]

A string representation of the generalized permutation.

INPUT:

  • sep - (default: ‘n’) a separator for the two intervals

OUTPUT:

string – the string that represents the permutation

EXAMPLES:

sage: from surface_dynamics import *

For permutations of iet:

sage: p = iet.Permutation('a b c','c b a')
sage: p.str()
'a b c\nc b a'
sage: p.str(sep=' | ')
'a b c | c b a'

The permutation can be rebuilt from the standard string:

sage: p == iet.Permutation(p.str())
True

For permutations of li:

sage: p = iet.GeneralizedPermutation('a b b','c c a')
sage: p.str()
'a b b\nc c a'
sage: p.str(sep=' | ')
'a b b | c c a'

sage: iet.GeneralizedPermutation([0,0], [1,2,3,2,1,3])
0 0
1 2 3 2 1 3
sage: iet.GeneralizedPermutation([0,1,2,1,0,2], [3,3])
0 1 2 1 0 2
3 3

Again, the generalized permutation can be rebuilt from the standard string:

sage: p == iet.GeneralizedPermutation(p.str())
True

With flips:

sage: p = iet.GeneralizedPermutation('a a','b b',flips='a')
sage: print(p.str())
-a -a
 b  b
 sage: print(p.str('/'))
 -a -a/ b  b

Alignment with alphabet of different sizes:

sage: p = iet.Permutation('aa b ccc d', 'd b ccc aa')
sage: print(p.str())
aa b ccc d
d b ccc aa
sage: print(p.str(align='left'))
aa b ccc d
d  b ccc aa
sage: print(p.str(align='right'))
aa b ccc  d
 d b ccc aa

sage: p = iet.GeneralizedPermutation('aa fff b ccc b fff d', 'eee d eee ccc aa')
sage: print(p.str())
aa fff b ccc b fff d
eee d eee ccc aa
sage: print(p.str(align='left'))
aa  fff b   ccc b  fff d
eee d   eee ccc aa
sage: print(p.str(align='right'))
 aa fff   b ccc  b fff d
eee   d eee ccc aa
symmetric()[source]

Returns the symmetric permutation.

The symmetric permutation is the composition of the top-bottom inversion and the left-right inversion (which are geometrically orientation reversing).

There are two other symmetries of permutation which are accessible via the methods Permutation.left_right_inverse() and Permutation.top_bottom_inverse().

OUTPUT: a permutation

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c','c b a',reduced=True)
sage: p.symmetric() == p
True
sage: p = iet.Permutation('a b c d','c a d b',reduced=True)
sage: q = p.symmetric()
sage: q
a b c d
b d a c
sage: q1 = p.tb_inverse().lr_inverse()
sage: q2 = p.lr_inverse().tb_inverse()
sage: q == q1 and q == q2
True

It works for any type of permutations:

sage: p = iet.GeneralizedPermutation('a b b','c c a',flips='ab')
sage: p
-a -b -b
 c  c -a
sage: p.symmetric()
-a  c  c
-b -b -a
tb_inverse(inplace=False)

Return the top-bottom inverse.

You can use also use the shorter .tb_inverse().

There are two other symmetries of permutation which are accessible via the methods Permutation.left_right_inverse() and Permutation.symmetric().

OUTPUT: a permutation

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b','b a')
sage: p.top_bottom_inverse()
b a
a b
sage: p = iet.Permutation('a b','b a',reduced=True)
sage: p.top_bottom_inverse() == p
True
sage: p = iet.Permutation('a b c d','c d a b')
sage: p.top_bottom_inverse()
c d a b
a b c d
top_bottom_inverse(inplace=False)[source]

Return the top-bottom inverse.

You can use also use the shorter .tb_inverse().

There are two other symmetries of permutation which are accessible via the methods Permutation.left_right_inverse() and Permutation.symmetric().

OUTPUT: a permutation

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b','b a')
sage: p.top_bottom_inverse()
b a
a b
sage: p = iet.Permutation('a b','b a',reduced=True)
sage: p.top_bottom_inverse() == p
True
sage: p = iet.Permutation('a b c d','c d a b')
sage: p.top_bottom_inverse()
c d a b
a b c d
vertical_inverse(inplace=False)

Return the left-right inverse.

The left-right inverse of a permutation, is the permutation obtained by reversing the order of the underlying ordering.

You can also use the shorter .lr_inverse()

There are two other symmetries of permutation which are accessible via the methods Permutation.top_bottom_inverse() and Permutation.symmetric().

OUTPUT: a permutation

EXAMPLES:

sage: from surface_dynamics import *

For labelled permutations:

sage: p = iet.Permutation('a b c','c a b')
sage: p.left_right_inverse()
c b a
b a c
sage: p = iet.Permutation('a b c d','c d a b')
sage: p.left_right_inverse()
d c b a
b a d c

for reduced permutations:

sage: p = iet.Permutation('a b c','c a b',reduced=True)
sage: p.left_right_inverse()
a b c
b c a
sage: p = iet.Permutation('a b c d','c d a b',reduced=True)
sage: p.left_right_inverse()
a b c d
c d a b

for labelled quadratic permutations:

sage: p = iet.GeneralizedPermutation('a a','b b c c')
sage: p.left_right_inverse()
a a
c c b b

for reduced quadratic permutations:

sage: p = iet.GeneralizedPermutation('a a','b b c c',reduced=True)
sage: p.left_right_inverse() == p
True
class surface_dynamics.interval_exchanges.template.PermutationIET(intervals=None, alphabet=None, reduced=False, flips=None)[source]

Bases: Permutation

has_rauzy_move(winner, side='right')[source]

Test if a Rauzy move can be performed on this permutation.

EXAMPLES:

sage: from surface_dynamics import *

for labelled permutations:

sage: p = iet.Permutation('a b c','a c b',reduced=False)
sage: p.has_rauzy_move(0,'right')
True
sage: p.has_rauzy_move(0,'left')
False
sage: p.has_rauzy_move(1,'right')
True
sage: p.has_rauzy_move(1,'left')
False

for reduced permutations:

sage: p = iet.Permutation('a b c','a c b',reduced=True)
sage: p.has_rauzy_move(0,'right')
True
sage: p.has_rauzy_move(0,'left')
False
sage: p.has_rauzy_move(1,'right')
True
sage: p.has_rauzy_move(1,'left')
False
is_irreducible(return_decomposition=False)[source]

Test irreducibility.

A permutation p = (p0,p1) is reducible if: set(p0[:i]) = set(p1[:i]) for an i < len(p0)

OUTPUT:

  • a boolean

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c', 'c b a')
sage: p.is_irreducible()
True

sage: p = iet.Permutation('a b c', 'b a c')
sage: p.is_irreducible()
False

sage: p = iet.Permutation('a b c', 'c b a', flips=['a'])
sage: p.is_irreducible()
True
twin(i, pos)[source]

Return the twin of the interval in the interval i at position pos.

EXAMPLES:

sage: from surface_dynamics import *
sage: p = iet.Permutation('a b c e d', 'e b d a c')
sage: p.twin(0,0)
(1, 3)
sage: p.twin(0,1)
(1, 1)

sage: twin_top = [p.twin(0,i) for i in range(p.length_top())]
sage: twin_bot = [p.twin(1,i) for i in range(p.length_bottom())]
sage: p.twin_list() == [twin_top, twin_bot]
True
twin_list()[source]

Returns the twin list of self.

The twin list is the involution without fixed point associated to that permutation seen as two lines of symbols. As the domain is two lines, the position are 2-tuples (i,j) where i specifies the line and j the position in the line.

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c','c b a')
sage: p.twin_list()[0]
[(1, 2), (1, 1), (1, 0)]
sage: p.twin_list()[1]
[(0, 2), (0, 1), (0, 0)]

We may check that it is actually an involution without fixed point:

sage: t = p.twin_list()
sage: all(t[i][j] != (i,j) for i in range(2) for j in range(len(t[i])))
True
sage: all(t[t[i][j][0]][t[i][j][1]] == (i,j) for i in range(2) for j in range(len(t[i])))
True
class surface_dynamics.interval_exchanges.template.PermutationLI(intervals=None, alphabet=None, reduced=False, flips=None)[source]

Bases: Permutation

erase_marked_points()[source]

Return a permutation without marked points.

This method is not implemented for generalized permutations.

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.GeneralizedPermutation('a a b','b c c')
sage: p.stratum()
Q_0(-1^4)
sage: p.erase_marked_points()
a a b
b c c
sage: p = iet.GeneralizedPermutation('a d d a b','b c c')
sage: p.stratum()
Q_0(0, -1^4)
sage: p.erase_marked_points()
Traceback (most recent call last):
...
NotImplementedError: not yet implemented! Do it!
genus()[source]

Returns the genus of any suspension of self.

The genus g can be deduced from the profile (see profile()) p=(p_1,ldots,p_k) of self by the formula: 4g-4 = sum_{i=1}^k (p_i - 2).

EXAMPLES:

sage: from surface_dynamics import *

sage: iet.GeneralizedPermutation('a a b','b c c').genus()
0
sage: iet.GeneralizedPermutation((0,1,2,1,3),(4,3,4,2,0)).genus()
2
has_rauzy_move(winner, side='right')[source]

Test of Rauzy movability (with an eventual specified choice of winner)

A quadratic (or generalized) permutation is rauzy_movable type depending on the possible length of the last interval. It’s dependent of the length equation.

INPUT:

  • winner - the integer ‘top’ or ‘bottom’

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.GeneralizedPermutation('a a','b b')
sage: p.has_rauzy_move('top','right')
False
sage: p.has_rauzy_move('top','left')
False
sage: p.has_rauzy_move('bottom','right')
False
sage: p.has_rauzy_move('bottom','left')
False
sage: p = iet.GeneralizedPermutation('a a b','b c c')
sage: p.has_rauzy_move('top','right')
True
sage: p.has_rauzy_move('bottom','right')
True
sage: p.has_rauzy_move('top','left')
True
sage: p.has_rauzy_move('bottom','left')
True
sage: p = iet.GeneralizedPermutation('a a','b b c c')
sage: p.has_rauzy_move('top','right')
True
sage: p.has_rauzy_move('bottom','right')
False
sage: p.has_rauzy_move('top','left')
True
sage: p.has_rauzy_move('bottom','left')
False
sage: p = iet.GeneralizedPermutation('a a b b','c c')
sage: p.has_rauzy_move('top','right')
False
sage: p.has_rauzy_move('bottom','right')
True
sage: p.has_rauzy_move('top','left')
False
sage: p.has_rauzy_move('bottom','left')
True
is_cylindric()[source]

Test if the permutation is cylindric

EXAMPLES:

sage: from surface_dynamics import *

sage: q = iet.GeneralizedPermutation('a b b','c c a')
sage: q.is_cylindric()
True
sage: q = iet.GeneralizedPermutation('a a b b','c c')
sage: q.is_cylindric()
False
is_hyperelliptic(verbose=False)[source]

Test if this permutation is in an hyperelliptic connected component.

EXAMPLES:

sage: from surface_dynamics import *

An example of hyperelliptic permutation:

sage: p = iet.GeneralizedPermutation([0,1,2,0,6,5,3,1,2,3],[4,5,6,4])
sage: p.is_hyperelliptic()
True

Check for the correspondence:

sage: q = Stratum([6,6], k=2)
sage: c_hyp, c_reg, c_irr = q.components()

sage: p_hyp = c_hyp.permutation_representative()
sage: p_hyp
0 1 2 3 4 1 5 6 7
7 6 5 8 4 3 2 8 0
sage: p_hyp.is_hyperelliptic()
True

sage: p_reg = c_reg.permutation_representative()
sage: p_reg
0 1 2 3 4 5 2 6 7 5
1 4 6 8 7 8 3 0
sage: p_reg.is_hyperelliptic()
False

sage: p_irr = c_irr.permutation_representative()
sage: p_irr
0 1 2 3 4 3 5 6 7
1 6 8 4 2 7 5 8 0
sage: p_irr.is_hyperelliptic()
False

sage: q = Stratum([3,3,2], k=2)
sage: c_hyp, c_non_hyp = q.components()
sage: p_hyp = c_hyp.permutation_representative()
sage: p_hyp.is_hyperelliptic()
True
sage: p_non_hyp = c_non_hyp.permutation_representative()
sage: p_non_hyp.is_hyperelliptic()
False
sage: q = Stratum([5,5,2], k=2)
sage: c_hyp, c_non_hyp = q.components()
sage: p_hyp = c_hyp.permutation_representative()
sage: p_hyp.is_hyperelliptic()
True
sage: p_non_hyp = c_non_hyp.permutation_representative()
sage: p_non_hyp.is_hyperelliptic()
False
sage: q = Stratum([3,3,1,1], k=2)
sage: c_hyp, c_non_hyp = q.components()
sage: p_hyp = c_hyp.permutation_representative()
sage: p_hyp.is_hyperelliptic()
True
sage: p_non_hyp = c_non_hyp.permutation_representative()
sage: p_non_hyp.is_hyperelliptic()
False
is_irreducible(return_decomposition=False)[source]

Test of reducibility

A quadratic (or generalized) permutation is reducible if there exists a decomposition

\[ \begin{align}\begin{aligned}A1 u B1 | ... | B1 u A2\\A1 u B2 | ... | B2 u A2\end{aligned}\end{align} \]

where no corners is empty, or exactly one corner is empty and it is on the left, or two and they are both on the right or on the left. The definition is due to [BoiLan09] where they prove that the property of being irreducible is stable under Rauzy induction.

INPUT:

  • return_decomposition - boolean (default: False) - if True, and the permutation is reducible, returns also the blocks A1 u B1, B1 u A2, A1 u B2 and B2 u A2 of a decomposition as above.

OUTPUT:

If return_decomposition is True, returns a 2-uple (test,decomposition) where test is the preceding test and decomposition is a 4-uple (A11,A12,A21,A22) where:

A11 = A1 u B1 A12 = B1 u A2 A21 = A1 u B2 A22 = B2 u A2

EXAMPLES:

sage: from surface_dynamics import *

sage: GP = iet.GeneralizedPermutation

sage: GP('a a','b b').is_irreducible()
False
sage: GP('a a b','b c c').is_irreducible()
True
sage: GP('1 2 3 4 5 1','5 6 6 4 3 2').is_irreducible()
True
marked_profile()[source]

Returns the marked profile of self.

The marked profile of a generalized permutation is an integer partition and some additional data associated to the angles of conical singularities in the suspension. The partition, called the profile, is the list of angles divided by 2pi (see profile()). The additional is called the marking and may be of two different types.

If the left endpoint and the right endpoint of the interval associated to the permutation coincides, then the marking is of type 1 and the additional data consists of a couple (m,a) such that m is the angle of the conical singularity and a is the angle between the outgoing separatrix associated to the left endpoint and the incoming separatrix associated to the right endpoint. A marking of type one is denoted m | a.

If the left endpoint and the right endpoint are two different conical singularities in the suspension, then the marking is of type 2 and the data consists in a couple (m_l,m_r) where m_l (resp. m_r) is the conical angle of the singularity at the left endpoint (resp. right endpoint). A marking of type two is denoted m_l circ m_r

EXAMPLES:

sage: from surface_dynamics import *

All possible markings for the profile [1, 1, 1, 1]:

sage: p = iet.GeneralizedPermutation('a a b','b c c')
sage: p.marked_profile()
1o1 [1, 1, 1, 1]
sage: p = iet.GeneralizedPermutation('a a','b b c c')
sage: p.marked_profile()
1|0 [1, 1, 1, 1]

All possible markings for the profile [4, 4]:

sage: p = iet.GeneralizedPermutation('0 1 2 1 3','3 4 0 4 2')
sage: p.marked_profile()
4o4 [4, 4]

sage: p = iet.GeneralizedPermutation('0 1 2 1 3','4 3 2 0 4')
sage: p.marked_profile()
4|0 [4, 4]

sage: p = iet.GeneralizedPermutation('0 1 0 2 3 2','4 3 4 1')
sage: p.marked_profile()
4|1 [4, 4]

sage: p = iet.GeneralizedPermutation('0 1 2 3 2','4 3 4 1 0')
sage: p.marked_profile()
4|2 [4, 4]

sage: p = iet.GeneralizedPermutation('0 1 0 1','2 3 2 4 3 4')
sage: p.marked_profile()
4|3 [4, 4]
marking()[source]

Return the marking induced by the two sides of the interval

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.GeneralizedPermutation('0 1 2 3 4 3 5 6 7','1 6 8 4 2 7 5 8 0')
sage: p.marking()
8|7

sage: p = iet.GeneralizedPermutation('0 1 2 3 4 3 5 6 7','1 6 8 4 2 7 8 0 5')
sage: p.marking()
8o8
orientation_cover()[source]

Return the orientation cover of this permutation.

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.GeneralizedPermutation('a a b', 'b c c')
sage: c = p.orientation_cover()
sage: c
Covering of degree 2 of the permutation:
a a b
b c c
sage: c.stratum()
H_1(0^4)

sage: C = Stratum([3,2,2,1], k=2).unique_component()
sage: p = C.permutation_representative()
sage: c = p.orientation_cover()
sage: c.stratum()
H_6(4, 2, 1^4)
profile()[source]

Returns the profile of self.

The profile of a generalized permutation is the list (d_1, ldots, d_k) where (d_1 pi, ldots, d_k pi) is the list of angles of any suspension of that generalized permutation.

See also marked_profile().

EXAMPLES:

sage: from surface_dynamics import *

sage: p1 = iet.GeneralizedPermutation('a a b','b c c')
sage: p1.profile()
[1, 1, 1, 1]
sage: all(p.profile() == [1, 1, 1, 1] for p in p1.rauzy_diagram())
True

sage: p2 = iet.GeneralizedPermutation('0 1 2 1 3','4 3 4 2 0')
sage: p2.profile()
[4, 4]
sage: all(p.profile() == [4,4] for p in p2.rauzy_diagram())
True

sage: p3 = iet.GeneralizedPermutation('0 1 2 3 3','2 1 4 4 0')
sage: p3.profile()
[3, 3, 1, 1]
sage: all(p.profile() == [3, 3, 1, 1] for p in p3.rauzy_diagram())
True
stratum_component()[source]

Return the connected component of stratum in which self belongs to.

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.GeneralizedPermutation('a b b','c c a')
sage: p.stratum_component()
Q_0(-1^4)^c

Test the exceptional strata in genus 3:

sage: Q = Stratum([9,-1], k=2)
sage: p = Q.regular_component().permutation_representative()
sage: p.stratum_component()
Q_3(9, -1)^reg
sage: p = Q.irregular_component().permutation_representative()
sage: p.stratum_component()
Q_3(9, -1)^irr

sage: Q = Stratum([6,3,-1], k=2)
sage: p = Q.regular_component().permutation_representative()
sage: p.stratum_component()
Q_3(6, 3, -1)^reg
sage: p = Q.irregular_component().permutation_representative()
sage: p.stratum_component()
Q_3(6, 3, -1)^irr

sage: Q = Stratum([3,3,3,-1], k=2)
sage: p = Q.regular_component().permutation_representative()
sage: p.stratum_component()
Q_3(3^3, -1)^reg
sage: p = Q.irregular_component().permutation_representative()
sage: p.stratum_component()
Q_3(3^3, -1)^irr

Test the exceptional strata in genus 4:

sage: Q = Stratum([12], k=2)
sage: p = Q.regular_component().permutation_representative()
sage: p.stratum_component()
Q_4(12)^reg
sage: p = Q.irregular_component().permutation_representative()
sage: p.stratum_component()
Q_4(12)^irr

sage: Q = Stratum([9,3], k=2)
sage: p = Q.regular_component().permutation_representative()
sage: p.stratum_component()
Q_4(9, 3)^reg
sage: p = Q.irregular_component().permutation_representative()
sage: p.stratum_component()
Q_4(9, 3)^irr

sage: Q = Stratum([6,6], k=2)
sage: p = Q.hyperelliptic_component().permutation_representative()
sage: p.stratum_component()
Q_4(6^2)^hyp
sage: p = Q.regular_component().permutation_representative()
sage: p.stratum_component()
Q_4(6^2)^reg
sage: p = Q.irregular_component().permutation_representative()
sage: p.stratum_component()
Q_4(6^2)^irr

sage: Q = Stratum([6,3,3], k=2)
sage: p = Q.regular_component().permutation_representative()
sage: p.stratum_component()
Q_4(6, 3^2)^reg
sage: p = Q.irregular_component().permutation_representative()
sage: p.stratum_component()
Q_4(6, 3^2)^irr

sage: Q = Stratum([3,3,3,3], k=2)
sage: p = Q.hyperelliptic_component().permutation_representative()
sage: p.stratum_component()
Q_4(3^4)^hyp
sage: p = Q.regular_component().permutation_representative()
sage: p.stratum_component()
Q_4(3^4)^reg
sage: p = Q.irregular_component().permutation_representative()
sage: p.stratum_component()
Q_4(3^4)^irr
to_cylindric()[source]

Return a cylindric permutation in the same extended Rauzy class

A generalized permutation is cylindric if the first letter in the top interval is the same as the last letter in the bottom interval or if the laster letter of the top interval is the same as the fist letter of the bottom interval.

EXAMPLES:

sage: from surface_dynamics import iet

sage: p = iet.GeneralizedPermutation('a b d a c','c e b e d')
sage: p.is_irreducible()
True
sage: p.to_cylindric().is_cylindric()
True

sage: p = iet.GeneralizedPermutation([0,1,2,1,2,3,4,5,6], [6,0,5,4,7,3,7], reduced=True)
sage: p.to_cylindric()
0 1 2 1 2 3 4 5 6
6 0 5 4 7 3 7

ALGORITHM:

The algorithm is naive. It computes the extended Rauzy class until it finds a cylindric permutation.

twin(i, pos)[source]

Return the twin of the letter in interval i at position pos

EXAMPLES:

sage: from surface_dynamics import *
sage: p = iet.GeneralizedPermutation('a a b c', 'c e b e')
sage: p.twin(0,0)
(0, 1)
sage: p.twin(0,1)
(0, 0)

sage: twin_top = [p.twin(0,i) for i in range(p.length_top())]
sage: twin_bot = [p.twin(1,i) for i in range(p.length_bottom())]
sage: p.twin_list() == [twin_top, twin_bot]
True
twin_list()[source]

Returns the twin list of self.

The twin list is the involution without fixed point which defines it. As the domain is naturally split into two lines we use a 2-tuple (i,j) to specify the element at position j in line i.

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.GeneralizedPermutation('a a b','b c c')
sage: p.twin_list()[0]
[(0, 1), (0, 0), (1, 0)]
sage: p.twin_list()[1]
[(0, 2), (1, 2), (1, 1)]

And we may check that it is actually an involution without fixed point:

sage: t = p.twin_list()
sage: all(t[i][j] != (i,j) for i in range(2) for j in range(len(t[i])))
True
sage: all(t[t[i][j][0]][t[i][j][1]] == (i,j) for i in range(2) for j in range(len(t[i])))
True

A slightly more complicated example:

sage: q = iet.GeneralizedPermutation('a b c a','d e f e g c b g d f')
sage: q.twin_list()[0]
[(0, 3), (1, 6), (1, 5), (0, 0)]
sage: q.twin_list()[1]
[(1, 8), (1, 3), (1, 9), (1, 1), (1, 7), (0, 2), (0, 1), (1, 4), (1, 0), (1, 2)]
sage: t = q.twin_list()
sage: all(t[t[i][j][0]][t[i][j][1]] == (i,j) for i in range(2) for j in range(len(t[i])))
True
class surface_dynamics.interval_exchanges.template.RauzyDiagram(p, right_induction=True, left_induction=False, left_right_inversion=False, top_bottom_inversion=False, symmetric=False)[source]

Bases: SageObject

Template for Rauzy diagrams.

AUTHORS:

  • Vincent Delecroix (2008-12-20): initial version

class Path(parent, *data)[source]

Bases: SageObject

Path in Rauzy diagram.

A path in a Rauzy diagram corresponds to a subsimplex of the simplex of lengths. This correspondence is obtained via the Rauzy induction. To a idoc IET we can associate a unique path in a Rauzy diagram. This establishes a correspondence between infinite full path in Rauzy diagram and equivalence topologic class of IET.

append(edge_type)[source]

Append an edge to the path.

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c','c b a')
sage: r = p.rauzy_diagram()
sage: g = r.path(p)
sage: g.append('top')
sage: g
Path of length 1 in a Rauzy diagram
sage: g.append('bottom')
sage: g
Path of length 2 in a Rauzy diagram
composition(function, composition=None)[source]

Compose an edges function on a path

INPUT:

  • path - either a Path or a tuple describing a path

  • function - function must be of the form

  • composition - the composition function

AUTHOR:

  • Vincent Delecroix (2009-09-29)

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b','b a')
sage: r = p.rauzy_diagram()
sage: def f(i,t):
....:     if t is None: return []
....:     return [t]
sage: g = r.path(p)
sage: g.composition(f,list.__add__)
[]
sage: g = r.path(p,0,1)
sage: g.composition(f, list.__add__)
[0, 1]
edge_types()[source]

Returns the edge types of the path.

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c','c b a')
sage: r = p.rauzy_diagram()
sage: g = r.path(p, 0, 1)
sage: g.edge_types()
[0, 1]
end()[source]

Returns the last vertex of the path.

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c','c b a')
sage: r = p.rauzy_diagram()
sage: g1 = r.path(p, 't', 'b', 't')
sage: g1.end() == p
True
sage: g2 = r.path(p, 'b', 't', 'b')
sage: g2.end() == p
True
extend(path)[source]

Extends self with another path.

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c d','d c b a')
sage: r = p.rauzy_diagram()
sage: g1 = r.path(p,'t','t')
sage: g2 = r.path(p.rauzy_move('t').rauzy_move('t'),'b','b')
sage: g = r.path(p,'t','t','b','b')
sage: g == g1 + g2
True
sage: g = copy(g1)
sage: g.extend(g2)
sage: g == g1 + g2
True
is_loop()[source]

Tests whether the path is a loop (start point = end point).

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b','b a')
sage: r = p.rauzy_diagram()
sage: r.path(p).is_loop()
True
sage: r.path(p,0,1,0,0).is_loop()
True
losers()[source]

Returns a list of the loosers on the path.

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c','c b a')
sage: r = p.rauzy_diagram()
sage: g0 = r.path(p,'t','b','t')
sage: g0.losers()
['a', 'c', 'b']
sage: g1 = r.path(p,'b','t','b')
sage: g1.losers()
['c', 'a', 'b']
pop()[source]

Pops the queue of the path

OUTPUT:

a path corresponding to the last edge

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b','b a')
sage: r = p.rauzy_diagram()
sage: g = r.path(p,0,1,0)
sage: g0,g1,g2,g3 = g[0], g[1], g[2], g[3]
sage: g.pop() == r.path(g2,0)
True
sage: g == r.path(g0,0,1)
True
sage: g.pop() == r.path(g1,1)
True
sage: g == r.path(g0,0)
True
sage: g.pop() == r.path(g0,0)
True
sage: g == r.path(g0)
True
sage: g.pop() == r.path(g0)
True
right_composition(function, composition=None)[source]

Compose an edges function on a path

INPUT:

  • function - function must be of the form (indice,type) -> element. Moreover function(None,None) must be an identity element for initialization.

  • composition - the composition function for the function. * if None (default: None)

start()[source]

Returns the first vertex of the path.

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c','c b a')
sage: r = p.rauzy_diagram()
sage: g = r.path(p, 't', 'b')
sage: g.start() == p
True
winners()[source]

Returns the winner list associated to the edge of the path.

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b','b a')
sage: r = p.rauzy_diagram()
sage: r.path(p).winners()
[]
sage: r.path(p,0).winners()
['b']
sage: r.path(p,1).winners()
['a']
alphabet(data=None)[source]
cardinality()[source]

Returns the number of permutations in this Rauzy diagram.

OUTPUT:

  • integer - the number of vertices in the diagram

EXAMPLES:

sage: from surface_dynamics import *

sage: r = iet.RauzyDiagram('a b','b a')
sage: r.cardinality()
1
sage: r = iet.RauzyDiagram('a b c','c b a')
sage: r.cardinality()
3
sage: r = iet.RauzyDiagram('a b c d','d c b a')
sage: r.cardinality()
7
complete(p)[source]

Completion of the Rauzy diagram.

Add to the Rauzy diagram all permutations that are obtained by successive operations defined by edge_types(). The permutation must be of the same type and the same length as the one used for the creation.

INPUT:

  • p - a permutation of Interval exchange transformation

Rauzy diagram is the reunion of all permutations that could be obtained with successive rauzy moves. This function just use the functions __getitem__ and has_rauzy_move and rauzy_move which must be defined for child and their corresponding permutation types.

edge_iterator()[source]

Returns an iterator over the edges of the graph.

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b','b a')
sage: r = p.rauzy_diagram()
sage: for e in r.edge_iterator():
....:     print('%s --> %s' %(e[0].str(sep='/'), e[1].str(sep='/')))
a b/b a --> a b/b a
a b/b a --> a b/b a
edge_to_loser(p=None, edge_type=None)[source]

Return the corresponding loser

edge_to_matrix(p=None, edge_type=None)[source]

Return the corresponding matrix

INPUT:

  • p - a permutation

  • edge_type - 0 or 1 corresponding to the type of the edge

OUTPUT:

A matrix

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c','c b a')
sage: d = p.rauzy_diagram()
sage: print(d.edge_to_matrix(p,1))
[1 0 1]
[0 1 0]
[0 0 1]
edge_to_winner(p=None, edge_type=None)[source]

Return the corresponding winner

edge_types()[source]

Print information about edges.

EXAMPLES:

sage: from surface_dynamics import *

sage: r = iet.RauzyDiagram('a b', 'b a')
sage: r.edge_types()
0: rauzy_move(0, -1)
1: rauzy_move(1, -1)
sage: r = iet.RauzyDiagram('a b', 'b a', left_induction=True)
sage: r.edge_types()
0: rauzy_move(0, -1)
1: rauzy_move(1, -1)
2: rauzy_move(0, 0)
3: rauzy_move(1, 0)
sage: r = iet.RauzyDiagram('a b',' b a',symmetric=True)
sage: r.edge_types()
0: rauzy_move(0, -1)
1: rauzy_move(1, -1)
2: symmetric()
edge_types_index(data)[source]

Try to convert the data as an edge type.

INPUT:

  • data - a string

OUTPUT:

integer

EXAMPLES:

sage: from surface_dynamics import *

For a standard Rauzy diagram (only right induction) the 0 index corresponds to the ‘top’ induction and the index 1 corresponds to the ‘bottom’ one:

sage: p = iet.Permutation('a b c','c b a')
sage: r = p.rauzy_diagram()
sage: r.edge_types_index('top')
0
sage: r[p][0] == p.rauzy_move('top')
True
sage: r.edge_types_index('bottom')
1
sage: r[p][1] == p.rauzy_move('bottom')
True

The special operations (inversion and symmetry) always appears after the different Rauzy inductions:

sage: p = iet.Permutation('a b c','c b a')
sage: r = p.rauzy_diagram(symmetric=True)
sage: r.edge_types_index('symmetric')
2
sage: r[p][2] == p.symmetric()
True

This function always try to resolve conflictuous name. If it’s impossible a ValueError is raised:

sage: p = iet.Permutation('a b c','c b a')
sage: r = p.rauzy_diagram(left_induction=True)
sage: r.edge_types_index('top')
Traceback (most recent call last):
...
ValueError: left and right inductions must be differentiated
sage: r.edge_types_index('top_right')
0
sage: r[p][0] == p.rauzy_move(0)
True
sage: r.edge_types_index('bottom_left')
3
sage: r[p][3] == p.rauzy_move('bottom', 'left')
True
sage: p = iet.Permutation('a b c','c b a')
sage: r = p.rauzy_diagram(left_right_inversion=True,top_bottom_inversion=True)
sage: r.edge_types_index('inversion')
Traceback (most recent call last):
...
ValueError: left-right and top-bottom inversions must be differentiated
sage: r.edge_types_index('lr_inverse')
2
sage: p.lr_inverse() == r[p][2]
True
sage: r.edge_types_index('tb_inverse')
3
sage: p.tb_inverse() == r[p][3]
True

Short names are accepted:

sage: p = iet.Permutation('a b c','c b a')
sage: r = p.rauzy_diagram(right_induction='top',top_bottom_inversion=True)
sage: r.edge_types_index('top_rauzy_move')
0
sage: r.edge_types_index('t')
0
sage: r.edge_types_index('tb')
1
sage: r.edge_types_index('inversion')
1
sage: r.edge_types_index('inverse')
1
sage: r.edge_types_index('i')
1
edges(labels=True)[source]

Returns a list of the edges.

EXAMPLES:

sage: from surface_dynamics import *

sage: r = iet.RauzyDiagram('a b','b a')
sage: len(r.edges())
2
graph()[source]

Returns the Rauzy diagram as a Graph object

The graph returned is more precisely a DiGraph (directed graph) with loops and multiedges allowed.

EXAMPLES:

sage: from surface_dynamics import *

sage: r = iet.RauzyDiagram('a b c','c b a')
sage: r
Rauzy diagram with 3 permutations
sage: r.graph()
Looped digraph on 3 vertices
letters()[source]

Returns the letters used by the RauzyDiagram.

EXAMPLES:

sage: from surface_dynamics import *

sage: r = iet.RauzyDiagram('a b','b a')
sage: r.alphabet()
{'a', 'b'}
sage: r.letters()
['a', 'b']
sage: r.alphabet('ABCDEF')
sage: r.alphabet()
{'A', 'B', 'C', 'D', 'E', 'F'}
sage: r.letters()
['A', 'B']
path(*data)[source]

Returns a path over this Rauzy diagram.

INPUT:

  • initial_vertex - the initial vertex (starting point of the path)

  • data - a sequence of edges

EXAMPLES:

sage: from surface_dynamics import *

sage: p = iet.Permutation('a b c','c b a')
sage: r = p.rauzy_diagram()
sage: g = r.path(p, 'top', 'bottom')
vertex_iterator()[source]

Returns an iterator over the vertices

EXAMPLES:

sage: from surface_dynamics import *

sage: r = iet.RauzyDiagram('a b','b a')
sage: for p in r.vertex_iterator(): print(p)
a b
b a
sage: r = iet.RauzyDiagram('a b c d','d c b a')
sage: r_1n = filter(lambda x: x.is_standard(), r)
sage: for p in r_1n: print(p)
a b c d
d c b a
vertices()[source]

Returns a list of the vertices.

EXAMPLES:

sage: from surface_dynamics import *

sage: r = iet.RauzyDiagram('a b','b a')
sage: for p in r.vertices(): print(p)
a b
b a
surface_dynamics.interval_exchanges.template.cylindric_canonical(p)[source]
surface_dynamics.interval_exchanges.template.interval_conversion(interval=None)[source]

Converts the argument in 0 or 1.

INPUT:

  • winner - ‘top’ (or ‘t’ or 0) or bottom (or ‘b’ or 1)

OUTPUT:

integer – 0 or 1

surface_dynamics.interval_exchanges.template.labelize_flip(couple)[source]

Returns a string from a 2-uple couple of the form (name, flip).

surface_dynamics.interval_exchanges.template.side_conversion(side=None)[source]

Converts the argument in 0 or -1.

INPUT:

  • side - either ‘left’ (or ‘l’ or 0) or ‘right’ (or ‘r’ or -1)

OUTPUT:

integer – 0 or -1

surface_dynamics.interval_exchanges.template.to_fat_graphs(twin)[source]