Triangulation

The main snappy class, namely Manifold, is derived from more basic class below.

class snappy.Triangulation

A Triangulation object represents a compact 3-manifold with torus boundary components, given as an ideal triangulation of the manifold’s interior. A Dehn-filling can be specified for each boundary component, allowing the description of closed 3-manifolds and some orbifolds. For non-orientable 3-manifolds, the boundary components can also be Klein bottles. Two Triangulations are equal (‘==’) if they represent combinatorially isomorphic triangulations. A Triangulation does not have any geometric structure, and usually one works with the subclass Manifold which adds this. Here’s a quick example:

>>> M = Triangulation('9_42')
>>> M.num_tetrahedra()
5
>>> M.is_orientable()
True

A Triangulation can be specified in a number of ways, e.g.

  • Triangulation(‘9_42’) : The complement of the knot 9_42 in S^3.

  • Triangulation(‘m125(1,2)(4,5)’) : The SnapPea census manifold m125 where the first cusp has Dehn filling (1,2) and the second cusp has filling (4,5).

  • Triangulation() : Opens a link editor window where can you specify a link complement.

In general, the specification can be from among the below, with information on Dehn fillings added.

  • SnapPea cusped census manifolds: e.g. ‘m123’, ‘s123’, ‘v123’.

  • Link complements:
    • Rolfsen’s table: e.g. ‘4_1’, ‘04_1’, ‘5^2_6’, ‘6_4^7’, ‘L20935’, ‘l104001’.

    • Knots and links up to 14 crossings from tabulations by Hoste and Thistlethwaite: e.g. ‘K12a456’ or ‘L13n579’.

    • Hoste-Thistlethwaite Knotscape table: e.g. ‘11a17’ or ‘12n345’

    • Dowker-Thistlethwaite code: e.g. ‘DT:[(6,8,2,4)]’, ‘DT:dadbcda’

  • Once-punctured torus bundles: e.g. ‘b++LLR’, ‘b+-llR’, ‘bo-RRL’, ‘bn+LRLR’

  • Fibered manifold associated to a braid: ‘Braid:[1,2,-3,4]’

    Here, the braid is thought of as a mapping class of the punctured disc, and this manifold is the corresponding mapping torus. If you want the braid closure, do (1,0) filling of the last cusp.

  • From mapping class group data using Twister:

    ‘Bundle(S_{1,1}, [a0, B1])’ or ‘Splitting(S_{1,0}, [b1, A0], [a0,B1])’

    See the help for the ‘twister’ module for more.

  • A SnapPea triangulation or link projection file: ‘filename’

    The file will be loaded if found in the current directory or the path given by the shell variable SNAPPEA_MANIFOLD_DIRECTORY.

  • A Regina-style isomorphism signature, such as ‘dLQbcccdxwb’.

  • A string containing the contents of a SnapPea triangulation or link projection file.

DT_code(alpha=False, flips=False)

Return the Dowker-Thistlethwaite code of this link complement, if it is a link complement. The DT code is intended to be an immutable attribute, for use with knot and link exteriors only, which is set only when the manifold was created.

Here is the Whitehead link:

>>> M = Manifold('L5a1')
>>> M.DT_code()
[(6, 8), (2, 10, 4)]
>>> M.DT_code(alpha=True)
'ebbccdaeb'
>>> M.DT_code(alpha=True, flips=True)
'ebbccdaeb.01110'
>>> M.DT_code(flips=True)
([(6, 8), (2, 10, 4)], [0, 1, 1, 1, 0])
alexander_polynomial(**kwargs)

Computes the multivariable Alexander polynomial of the manifold:

sage: M = Manifold('K12n123')
sage: M.alexander_polynomial()
2*a^6 - 14*a^5 + 34*a^4 - 45*a^3 + 34*a^2 - 14*a + 2

sage: N = Triangulation('v1539(5,1)')
sage: N.alexander_polynomial()
a^2*b + a*b^2 + a*b + a + b

Any provided keyword arguments are passed to fundamental_group and so affect the group presentation used in the computation.

copy()

Returns a copy of the triangulation.

>>> M = Triangulation('m125')
>>> N = M.copy()
cover(permutation_rep)

Returns a Triangulation representing the finite cover specified by a transitive permutation representation. The representation is specified by a list of permutations, one for each generator of the simplified presentation of the fundamental group. Each permutation is specified as a list P such that set(P) == set(range(d)) where d is the degree of the cover.

>>> M = Triangulation('m004')
>>> N0 = M.cover([[1, 3, 0, 4, 2], [0, 2, 1, 4, 3]])
>>> N0.homology()
Z + Z + Z
>>> N0.cover_info()['type']
'irregular'
>>> N0.cover_info()['base']
'm004'
>>> N0.cover_info()['degree']
5

Within Sage the permutations can also be of type PermutationGroupElement, in which case they act on the set range(1, d + 1). Or, you can specify a GAP or Magma subgroup of the fundamental group. For examples, see the docstring for Manifold.cover

cover_info()

If this is a manifold or triangulation which was constructed as a covering space, return a dictionary describing the cover. Otherwise return 0. The dictionary keys are ‘base’, ‘type’ and ‘degree’.

covers(degree, method=None, cover_type='all')

M.covers(degree, method=None, cover_type=’all’)

Returns a list of Triangulations corresponding to all of the finite covers of the given degree. The default method is ‘low_index’ for general covers and ‘snappea’ for cyclic covers. The former uses Sim’s algorithm while the latter uses the original Snappea algorithm.

WARNING: If the degree is large this might take a very, very, very long time.

>>> M = Triangulation('m003')
>>> covers = M.covers(4)
>>> sorted(N.homology() for N in covers)
[Z/3 + Z/15 + Z, Z/5 + Z + Z]

It is faster to look just at cyclic covers.

>>> covers = M.covers(4, cover_type='cyclic')
>>> [(N, N.homology()) for N in covers]
[(m003~cyc~0(0,0), Z/3 + Z/15 + Z)]

Here we check that we get the same number of covers with the ‘snappea’ and ‘low_index’ methods.

>>> M = Triangulation('m125')
>>> len(M.covers(5))
19
>>> len(M.covers(5, method='snappea'))
19

If you are using Sage, you can use GAP to find the subgroups, which is often much faster, by specifying the optional argument method = ‘gap’ If you have Magma installed, you can used it to do the heavy lifting by specifying method=’magma’.

cusp_info(data_spec=None)

Returns an info object containing information about the given cusp. Usage:

>>> M = Triangulation('v3227(0,0)(1,2)(3,2)')
>>> M.cusp_info(1)
Cusp 1 : torus cusp with Dehn filling coefficients (M, L) = (1.0, 2.0)
>>> c = M.cusp_info(1)
>>> c.is_complete
False
>>> sorted(c.keys())
['filling', 'index', 'is_complete', 'topology']

You can get information about multiple cusps at once:

>>> M.cusp_info()
[Cusp 0 : torus cusp, not filled,
 Cusp 1 : torus cusp with Dehn filling coefficients (M, L) = (1.0, 2.0),
 Cusp 2 : torus cusp with Dehn filling coefficients (M, L) = (3.0, 2.0)]
>>> M.cusp_info('is_complete')
[True, False, False]
dehn_fill(filling_data, which_cusp=None)

Set the Dehn filling coefficients of the cusps. This can be specified in the following ways, where the cusps are numbered by 0,1,…,(num_cusps - 1).

  • Fill cusp 2:

    >>> M = Triangulation('8^4_1')
    >>> M.dehn_fill((2,3), 2)
    >>> M
    8^4_1(0,0)(0,0)(2,3)(0,0)
    
  • Fill the last cusp:

    >>> M.dehn_fill((1,5), -1)
    >>> M
    8^4_1(0,0)(0,0)(2,3)(1,5)
    
  • Fill the first two cusps:

    >>> M.dehn_fill( [ (3,0), (1, -4) ])
    >>> M
    8^4_1(3,0)(1,-4)(2,3)(1,5)
    
  • When there is only one cusp, there’s a shortcut

    >>> N = Triangulation('m004')
    >>> N.dehn_fill( (-3,4) )
    >>> N
    m004(-3,4)
    

Does not return a new Triangulation.

edge_valences()

Returns a dictionary whose keys are the valences of the edges in the triangulation, and the value associated to a key is the number of edges of that valence.

>>> M = Triangulation('v3227')
>>> M.edge_valences()     
{10: 1, 4: 1, 5: 2, 6: 3}

For a triangulation of the exterior of a link in the 3-sphere, return a planar diagram for the link. The peripheral curves whose Dehn filling is the 3-sphere are part of the input, specified by either:

  1. If no cusp is filled, then they are the meridians of the current peripheral curves.

  2. If every cusp is filled, then they are the current Dehn filling curves.

In particular, it does not try to determine whether there exist fillings on the input which give the 3-sphere. Example usage:

>>> M = Manifold('m016')
>>> L = exterior_to_link(M)
>>> L.exterior().is_isometric_to(M)
True

The algorithm used is that of Dunfield, Obeidin, and Rudd. The optional arguments are as follows.

  • verbose: When True, prints progress updates as the algorithm goes along.

  • check_input: When True (the default), first checks that the fundamental group of the specified Dehn filling is trivial. As it doesn’t try too hard to simplify the group presentation, it can happen that this check fails but the algorithm still finds a diagram if you pass check_input=False.

  • check_answer: When True (the default), take the exterior of the final link diagram and use Manifold.is_isometric_to to confirm that it is homeomorphic to the input. If the input is not hyperbolic or is very large, this check may fail even though the diagram is correct.

  • careful_perturbation: The rational coordinates of the intermediate PL links are periodically rounded to control the size of their denominators. When careful_perturbation=True (the default), computations are performed to ensure this rounding does not change the isotopy class of the link.

  • simplify_link: When True (the default), uses Link.simplify('global') to minimize the size of the final diagram; otherwise, it just does basic simplifications, which can be much faster if the initial link is complicated.

  • pachner_search_tries: Controls how hard to search for a suitable sequence of Pachner moves from the filled input triangulation to a standard triangulation of the 3-sphere.

  • seed: The algorithm involves many random choices, and hence each run typically produces a different diagram of the underlying link. If you need the same output each time, you can specify a fixed seed for the various pseudo-random number generators.

Note on rigor: Provided at least one of check_answer and careful_perturbation is True, the exterior of the output link is guaranteed to match the input (including the choice of meridians).

Warning: The order of the link components and the cusps of the input manifold is only guaranteed to match when check_answer=True. Even then, the implicit orientation along each component of the link may not be preserved.

filled_triangulation(cusps_to_fill='all')

Return a new manifold where the specified cusps have been permanently filled in. Examples:

Filling all the cusps:

>>> M = Triangulation('m125(1,2)(3,4)')
>>> N = M.filled_triangulation()
>>> N.num_cusps()
0

Filling cusps 0 and 2 :

>>> M = Triangulation('v3227(1,2)(3,4)(5,6)')
>>> M.filled_triangulation([0,2])
v3227_filled(3,4)
fundamental_group(simplify_presentation=True, fillings_may_affect_generators=True, minimize_number_of_generators=True, try_hard_to_shorten_relators=True)

Returns a FundamentalGroup object representing the fundamental group of the manifold. If integer Dehn surgery parameters have been set, then the corresponding peripheral elements are killed.

>>> M = Triangulation('m004')
>>> G = M.fundamental_group()
>>> G
Generators:
   a,b
Relators:
   aaabABBAb
>>> G.peripheral_curves()
[('ab', 'aBAbABab')]

There are four optional arguments all of which default to True:

  • simplify_presentation

  • fillings_may_affect_generators

  • minimize_number_of_generators

  • try_hard_to_shorten_relators

>>> M.fundamental_group(False, False, False)
Generators:
   a,b,c
Relators:
   CbAcB
   BacA
gluing_equations(form='log')

In the default mode, this function returns a matrix with rows of the form

a b c d e f …

which means

a*log(z0) + b*log(1/(1-z0)) + c*log((z0-1)/z0) + d*log(z1) +… = 2 pi i

for an edge equation, and (same) = 0 for a cusp equation. Here, the cusp equations come at the bottom of the matrix, and are listed in the form: meridian of cusp 0, longitude of cusp 0, meridian of cusp 1, longitude of cusp 1,…

In terms of the tetrahedra, a is the invariant of the edge (2,3), b the invariant of the edge (0,2) and c is the invariant of the edge (1,2). See kernel_code/edge_classes.c for a detailed account of the convention used.

If the optional argument form=’rect’ is given, then this function returns a list of tuples of the form:

( [a0, a1,..,a_n], [b_0, b_1,…,b_n], c)

where this corresponds to the equation

z0^a0 (1 - z0)^b0 z1^a1(1 - z1)^b1 … = c

where c = 1 or -1.

>>> M = Triangulation('m004(2,3)')
>>> M.gluing_equations()
[ 2  1  0  1  0  2]
[ 0  1  2  1  2  0]
[ 2  0  0  0 -8  6]
>>> M.gluing_equations(form='rect')
[([2, -1], [-1, 2], 1), ([-2, 1], [1, -2], 1), ([2, -6], [0, 14], 1)]
gluing_equations_pgl(N=2, equation_type='all')

M.gluing_equations_pgl(N = 2, equation_type=’all’)

Returns a NeumannZagierTypeEquations object that contains a matrix encoding the gluing equations for boundary-parabolic PGL(N,C) representations together with explanations of the meaning of the rows and the columns of the matrix.

This method generalizes gluing_equations() to PGL(N,C)-representations as described in Stavros Garoufalidis, Matthias Goerner, Christian K. Zickert: “Gluing Equations for PGL(n,C)-Representations of 3-Manifolds” (http://arxiv.org/abs/1207.6711).

The result of the traditional gluing_equations() can be obtained from the general method by:

>>> M = Triangulation('m004')
>>> M.gluing_equations_pgl().matrix
[ 2  1  0  1  0  2]
[ 0  1  2  1  2  0]
[ 1  0  0  0 -1  0]
[ 0  0  0  0 -2  2]

But besides the matrix, the method also returns explanations of the columns and rows:

>>> M = Triangulation("m004")
>>> M.gluing_equations_pgl()
NeumannZagierTypeEquations(
  [ 2  1  0  1  0  2]
  [ 0  1  2  1  2  0]
  [ 1  0  0  0 -1  0]
  [ 0  0  0  0 -2  2],
  explain_columns = ['z_0000_0', 'zp_0000_0', 'zpp_0000_0', 'z_0000_1', 'zp_0000_1', 'zpp_0000_1'],
  explain_rows = ['edge_0_0', 'edge_0_1', 'meridian_0_0', 'longitude_0_0'])

The first row of the matrix means that the edge equation for edge 0 is

z_0000_0 ^ 2 * zp_0000_0 * z_0000_1 * zpp_0000_1 ^ 2 = 1.

Similarly, the next row encodes the edge equation for the other edge and the next two rows encode peripheral equations.

Following the SnapPy convention, a z denotes the cross ratio z at the edge (0,1), a zp the cross ratio z’ at the edge (0,2) and a zpp the cross ratio z” at the edge (1,2). The entire symbol z_xxxx_y then denotes the cross ratio belonging to the subsimplex at integral point xxxx (always 0000 for N = 2) of the simplex y. Note: the SnapPy convention is different from the paper mentioned above, e.g., compare kernel_code/edge_classes.c with Figure 3. We follow the SnapPy convention here so that all computations done in SnapPy are consistent.

The explanations of the rows and columns can be obtained explicitly by:

>>> M.gluing_equations_pgl(N = 3, equation_type = 'peripheral').explain_rows
['meridian_0_0', 'meridian_1_0', 'longitude_0_0', 'longitude_1_0']
>>> M.gluing_equations_pgl(N = 2).explain_columns
['z_0000_0', 'zp_0000_0', 'zpp_0000_0', 'z_0000_1', 'zp_0000_1', 'zpp_0000_1']

A subset of all gluing equations can be obtained by setting the equation_type:

  • all gluing equations: ‘all’

  • non-peripheral equations: ‘non_peripheral’

    • edge gluing equations: ‘edge’

    • face gluing equations: ‘face’

    • internal gluing equations: ‘internal’

  • cusp gluing equations: ‘peripheral’

    • cusp gluing equations for meridians: ‘meridian’

    • cusp gluing equations for longitudes: ‘longitude’

has_finite_vertices()

Returns True if and only if the triangulation has finite (non-ideal) vertices.

>>> T = Triangulation("m004")
>>> T.has_finite_vertices()
False
>>> T.dehn_fill((12,13))
>>> S = T.filled_triangulation()
>>> S.has_finite_vertices()
True

When trying to find a hyperbolic structure, SnapPea will eliminate finite vertices:

>>> M = S.with_hyperbolic_structure()
>>> M.has_finite_vertices()
False
homological_longitude(cusp=None)

Returns the peripheral curve in the given cusp, if any, which is homologically trivial (with rational coefficients) in the manifold:

sage: M = Manifold('m015')
sage: M.homological_longitude()
(2, -1)

If no cusp is specified, the default is the first unfilled cusp; if all cusps are filled, the default is the first cusp:

sage: M = Manifold('L5a1(3,4)(0,0)')
sage: M.homological_longitude()
(0, 1)

The components of the next link have nontrivial linking number so there is no such curve:

sage: W = Manifold('L7a2')
sage: W.homological_longitude(cusp=1) is None
True

If every curve in the given cusp is trivial in the rational homology of the manifold, an exception is raised:

sage: M = Manifold('4_1(1,0)')
sage: M.homological_longitude()
Traceback (most recent call last):
...
ValueError: Every curve on cusp is homologically trivial
homology()

Returns an AbelianGroup representing the first integral homology group of the underlying (Dehn filled) manifold.

>>> M = Triangulation('m003')
>>> M.homology()
Z/5 + Z
is_orientable()

Return whether the underlying 3-manifold is orientable.

>>> M = Triangulation('x124')
>>> M.is_orientable()
False
isomorphisms_to(Triangulation other)

Returns a complete list of combinatorial isomorphisms between the two triangulations:

>>> M = Manifold('5^2_1')
>>> N = Manifold('5^2_1')
>>> N.set_peripheral_curves([[[2,3],[-1,-1]],[[1,1],[0,1]]])
>>> isoms = M.isomorphisms_to(N)
>>> isoms[6]
0 -> 1  1 -> 0
[ 1 0]  [-1 1]
[-1 1]  [-3 2]
Does not extend to link

Each transformation between cusps is given by a matrix which acts on the left. That is, the two columns of the matrix give the image of the meridian and longitude respectively. In the above example, the meridian of cusp 0 is sent to the meridian of cusp 1.

If the manifold is stored as a link complement in your current session then it returns the number of components and crossing of the link. To view and interact with the link see spherogram.Link.view() and Manifold.plink().

name()

Return the name of the triangulation.

>>> M = Triangulation('4_1')
>>> M.name()
'4_1'
normal_boundary_slopes(subset='all', algorithm='FXrays')

For a one-cusped manifold, returns all the nonempty boundary slopes of spun normal surfaces. Provided the triangulation supports a genuine hyperbolic structure, then by Thurston and Walsh any strict boundary slope (the boundary of an essential surface which is not a fiber or semifiber) must be listed here.

>>> M = Manifold('K3_1')
>>> M.normal_boundary_slopes()
[(16, -1), (20, -1), (37, -2)]

If the subset flag is set to 'kabaya', then it only returns boundary slopes associated to vertex surfaces with a quad in every tetrahedron; by Theorem 1.1. of [DG] these are all strict boundary slopes.

>>> N = Manifold('m113')
>>> N.normal_boundary_slopes()
[(1, 1), (1, 2), (2, -1), (2, 3), (8, 11)]
>>> N.normal_boundary_slopes('kabaya')
[(8, 11)]

If the subset flag is set to 'brasile' then it returns only the boundary slopes that are associated to vertex surfaces giving isolated rays in the space of embedded normal surfaces.

>>> N.normal_boundary_slopes('brasile')
[(1, 2), (8, 11)]
normal_surfaces(algorithm='FXrays')

All the vertex spun-normal surfaces in the current triangulation.

>>> M = Manifold('m004')
>>> M.normal_surfaces()    
[<Surface 0: [0, 0] [1, 2] (4, 1)>,
 <Surface 1: [0, 1] [1, 2] (4, -1)>,
 <Surface 2: [1, 2] [2, 1] (-4, -1)>,
 <Surface 3: [2, 2] [2, 1] (-4, 1)>]
num_cusps(cusp_type='all')

Return the total number of cusps. By giving the optional argument ‘orientable’ or ‘nonorientable’ it will only count cusps of that type.

>>> M = Triangulation('m125')
>>> M.num_cusps()
2
num_tetrahedra()

Return the number of tetrahedra in the triangulation.

>>> M = Triangulation('m004')
>>> M.num_tetrahedra()
2
orientation_cover()

For a non-orientable Triangulation, returns the 2-fold cover which is orientable.

>>> X = Triangulation('x123')
>>> Y = X.orientation_cover()
>>> (X.is_orientable(), Y.is_orientable())
(False, True)
>>> Y
x123~(0,0)(0,0)
>>> Y.cover_info()['type']
'cyclic'
pickle()

Brings up a link editor window if the manifold is stored as a link complement in your current session.

>>> M = Manifold('4_1') # stored as a triangulation with a link
>>> M.link()
<Link: 1 comp; 4 cross>
>>> N = Manifold('m004') # stored as a triangulation without a link
>>> N.link() 
Traceback (most recent call last):
...
ValueError: No associated link known.
ptolemy_generalized_obstruction_classes(N)

M.ptolemy_generalized_obstruction_classes(N)

Returns the obstruction classes needed to compute PGL(N,C)-representations for any N, i.e., it returns a list with a representative cocycle for each element in H^2(M, boundary M; Z/N) / (Z/N)^* where (Z/N)^* are the units in Z/N. The first element in the list always corresponds to the trivial obstruction class. The generalized ptolemy obstruction classes are thus a generalization of the ptolemy obstruction classes that allow to find all boundary-unipotent PGL(N,C)-representations including those that do not lift to boundary-unipotent SL(N,C)-representations for N odd or SL(N,C)/{+1,-1}-representations for N even.

For example, 4_1 has three obstruction classes up to equivalence:

>>> M = Manifold("4_1")
>>> c = M.ptolemy_generalized_obstruction_classes(4)
>>> len(c)
3

For 4_1, we only get three obstruction classes even though we have H^2(M, boundary M; Z/4) = Z/4 because the two obstruction classes 1 in Z/4 and -1 in Z/4 are related by a unit and thus give isomorphic Ptolemy varieties.

The primary use of an obstruction class sigma is to construct the Ptolemy variety of sigma. This variety computes boundary-unipotent PGL(N,C)-representations whose obstruction class to a boundary-unipotent lift to SL(N,C) is sigma.

For example for 4_1, there are 2 obstruction classes for N = 3:

>>> M = Manifold("4_1")
>>> c = M.ptolemy_generalized_obstruction_classes(3)
>>> len(c)
2

The Ptolemy variety parametrizing boundary-unipotent SL(3,C)-representations of 4_1 is obtained by

>>> p = M.ptolemy_variety(N = 3, obstruction_class = c[0])

and the Ptolemy variety parametrizing boundary-unipotent PSL(3,C)-representations of 4_1 that do not lift to boundary-unipotent SL(3,C)-representations is obtained by

>>> p = M.ptolemy_variety(N = 3, obstruction_class = c[1])

The cocycle representing the non-trivial obstruction class looks as follows:

>>> c[1]
PtolemyGeneralizedObstructionClass([2, 0, 0, 1])

This means that the cocycle takes the value -1 in Z/3 on the first face class and 1 on the fourth face class but zero on every other of the four face classes.

ptolemy_obstruction_classes()

Returns the obstruction classes needed to compute pSL(N,C) = SL(N,C)/{+1,-1} representations for even N, i.e., it returns a list with a representative cocycle for each class in H^2(M, boundary M; Z/2). The first element in the list is always representing the trivial obstruction class.

For example, 4_1 has two obstruction classes:

>>> M = Manifold("4_1")
>>> c = M.ptolemy_obstruction_classes()
>>> len(c)
2

The primary use of these obstruction classes is to construct the Ptolemy variety as described in Definition 1.7 of Stavros Garoufalidis, Dylan Thurston, Christian K. Zickert: “The Complex Volume of SL(n,C)-Representations of 3-Manifolds” (http://arxiv.org/abs/1111.2828).

For example, to construct the Ptolemy variety for PSL(2,C)-representations of 4_1 that do not lift to boundary-parabolic SL(2,C)-representations, use:

>>> p = M.ptolemy_variety(N = 2, obstruction_class = c[1])

Or the following short-cut:

>>> p = M.ptolemy_variety(2, obstruction_class = 1)

Note that this obstruction class only makes sense for even N:

>>> p = M.ptolemy_variety(3, obstruction_class = c[1])
Traceback (most recent call last):
...
AssertionError: PtolemyObstructionClass only makes sense for even N, try PtolemyGeneralizedObstructionClass

To obtain PGL(N,C)-representations for N > 2, use the generalized obstruction class:

>>> c = M.ptolemy_generalized_obstruction_classes(3)
>>> p = M.ptolemy_variety(3, obstruction_class = c[1])

The original obstruction class encodes a representing cocycle in Z/2 as follows:

>>> c = M.ptolemy_obstruction_classes()
>>> c[1]
PtolemyObstructionClass(s_0_0 + 1, s_1_0 - 1, s_2_0 - 1, s_3_0 + 1, s_0_0 - s_0_1, s_1_0 - s_3_1, s_2_0 - s_2_1, s_3_0 - s_1_1)

This means that the cocycle to represent this obstruction class in Z/2 takes value 1 in Z/2 on face 0 of tetrahedra 0 (because s_0_0 = -1) and value 0 in Z/2 on face 1 of tetrahedra 0 (because s_1_0 = +1).

Face 3 of tetrahedra 0 and face 1 of tetrahedra 1 are identified, hence the cocycle takes the same value on those two faces (s_3_0 = s_1_1).

ptolemy_variety(N, obstruction_class=None, simplify=True, eliminate_fixed_ptolemys=False)

M.ptolemy_variety(N, obstruction_class = None, simplify = True, eliminate_fixed_ptolemys = False)

Returns a Ptolemy variety as described in

  • Stavros Garoufalidis, Dyland Thurston, Christian K. Zickert: “The Complex Volume of SL(n,C)-Representations of 3-Manifolds” (http://arxiv.org/abs/1111.2828)

  • Stavros Garoufalidis, Matthias Goerner, Christian K. Zickert: “Gluing Equations for PGL(n,C)-Representations of 3-Manifolds ” (http://arxiv.org/abs/1207.6711)

The variety can be exported to magma or sage and solved there. The solutions can be processed to compute invariants. The method can also be used to automatically look up precomputed solutions from the database at http://ptolemy.unhyperbolic.org/data .

Example for m011 and PSL(2,C)-representations:

>>> M = Manifold("m011")

Obtain all Ptolemy varieties for PSL(2,C)-representations:

>>> p = M.ptolemy_variety(2, obstruction_class = 'all')

There are two Ptolemy varieties for the two obstruction classes:

>>> len(p)
2

Retrieve the solutions from the database

>>> sols = p.retrieve_solutions() 

Compute the solutions using magma (default in SnapPy)

>>> sols = p.compute_solutions(engine = 'magma') 

Compute the solutions using singular (default in sage)

>>> sols = p.compute_solutions(engine = 'sage') 

Note that magma is significantly faster.

Compute all resulting complex volumes

>>> cvols = sols.complex_volume_numerical() 
>>> cvols  
[[[-4.29405713186238 E-16 + 0.725471193740844*I,
   -0.942707362776931 + 0.459731436553693*I,
   0.942707362776931 + 0.459731436553693*I]],
 [[3.94159248086745 E-15 + 0.312682687518267*I,
   4.64549527022581 E-15 + 0.680993020093457*I,
   -2.78183391239608 - 0.496837853805869*I,
   2.78183391239608 - 0.496837853805869*I]]]

Show complex volumes as a non-nested list:

>>> cvols.flatten(depth=2) 
[-4.29405713186238 E-16 + 0.725471193740844*I,
 -0.942707362776931 + 0.459731436553693*I,
 0.942707362776931 + 0.459731436553693*I,
 3.94159248086745 E-15 + 0.312682687518267*I,
 4.64549527022581 E-15 + 0.680993020093457*I,
 -2.78183391239608 - 0.496837853805869*I,
 2.78183391239608 - 0.496837853805869*I]

For more examples, go to http://ptolemy.unhyperbolic.org/

=== Optional Arguments ===

obstruction_class — class from Definition 1.7 of (1). None for trivial class or a value returned from ptolemy_obstruction_classes. Short cuts: obstruction_class = ‘all’ returns a list of Ptolemy varieties for each obstruction. For easier iteration, can set obstruction_class to an integer.

simplify — boolean to indicate whether to simplify the equations which significantly reduces the number of variables. Simplifying means that several identified Ptolemy coordinates x = y = z = … are eliminated instead of adding relations x - y = 0, y - z = 0, …

eliminate_fixed_ptolemys — boolean to indicate whether to eliminate the Ptolemy coordinates that are set to 1 for fixing the decoration. Even though this simplifies the resulting representation, setting it to True can cause magma to run longer when finding a Groebner basis.

=== Examples for 4_1 ===

>>> M = Manifold("4_1")

Get the varieties for all obstruction classes at once (use help(varieties[0]) for more information):

>>> varieties = M.ptolemy_variety(2, obstruction_class = "all")

Print the variety as an ideal (sage object) for the non-trivial class:

>>> varieties[1].ideal    
Ideal (-c_0011_0^2 + c_0011_0*c_0101_0 + c_0101_0^2, -c_0011_0^2 - c_0011_0*c_0101_0 + c_0101_0^2, c_0011_0 - 1) of Multivariate Polynomial Ring in c_0011_0, c_0101_0 over Rational Field

Print the equations of the variety for the non-trivial class:

>>> for eqn in varieties[1].equations:
...     print(eqn)          
     - c_0011_0 * c_0101_0 + c_0011_0^2 + c_0101_0^2
     c_0011_0 * c_0101_0 - c_0011_0^2 - c_0101_0^2
     - 1 + c_0011_0

Generate a magma file to compute Primary Decomposition for N = 3:

>>> p = M.ptolemy_variety(3)
>>> s = p.to_magma()
>>> print(s.split("ring and ideal")[1].strip())     
R<c_0012_0, c_0012_1, c_0102_0, c_0111_0, c_0201_0, c_1011_0, c_1011_1, c_1101_0> := PolynomialRing(RationalField(), 8, "grevlex");
MyIdeal := ideal<R |
          c_0012_0 * c_1101_0 + c_0102_0 * c_0111_0 - c_0102_0 * c_1011_0,
    ...

=== If you have a magma installation ===

Call p.compute_solutions() to automatically call magma on the above output and produce exact solutions!!!

>>> try:
...     sols = p.compute_solutions()
... except:
...     # magma failed, use precomputed_solutions
...     sols = None

Check solutions against manifold >>> if sols: … dummy = sols.check_against_manifold()

=== If you do not have a magma installation ===

Load a precomputed example from magma which is provided with the package:

>>> from snappy.ptolemy.processMagmaFile import _magma_output_for_4_1__sl3, solutions_from_magma
>>> print(_magma_output_for_4_1__sl3)      

==TRIANGULATION=BEGINS==
% Triangulation
4_1
...

Parse the file and produce solutions:

>>> sols = solutions_from_magma(_magma_output_for_4_1__sl3)
>>> dummy = sols.check_against_manifold()

=== Continue here whether you have or do not have magma ===

Pick the first solution of the three different solutions (up to Galois conjugates):

>>> len(sols)
3
>>> solution = sols[0]

Read the exact value for c_1020_0 (help(solution) for more information on how to compute cross ratios, volumes and other invariants):

>>> solution['c_1020_0']
Mod(-1/2*x - 3/2, x^2 + 3*x + 4)

Example of simplified vs non-simplified variety for N = 4:

>>> simplified = M.ptolemy_variety(4, obstruction_class = 1)
>>> full = M.ptolemy_variety(4, obstruction_class = 1, simplify = False)
>>> len(simplified.variables), len(full.variables)
(21, 63)
>>> len(simplified.equations), len(full.equations)
(24, 72)
randomize(blowup_multiple=4, passes_at_fours=6)

Perform random Pachner moves on the underlying triangulation, including some initial 3 -> 2 moves that increase the number of tetrahedra by blowup_multiple.

>>> M = Triangulation('Braid:[1,2,-3,-3,1,2]')
>>> M.randomize()
reverse_orientation()

Reverses the orientation of the Triangulation, presuming that it is orientable.

>>> M = Manifold('m015')
>>> cs = M.chern_simons()
>>> M.reverse_orientation()
>>> abs(cs + M.chern_simons()) 
0.0
save(file_name=None)

Save the triangulation as a SnapPea triangulation file.

>>> M = Triangulation('m004')
>>> M.save('fig-eight.tri')     

To retrieve a SnapPea triangulation from the saved file you can do the following. The first command creates a cusped manifold M. The second one creates the filled manifold M1 with Dehn coefficients (2,3).

>>> M = Manifold('fig-eight.tri')   
>>> M1 = Manifold('fig-eight.tri(2,3)')   
set_name(new_name)

Give the triangulation a new name.

>>> M = Triangulation('4_1')
>>> M.set_name('figure-eight-comp')
>>> M
figure-eight-comp(0,0)
set_peripheral_curves(peripheral_data, which_cusp=None, return_matrices=False)

Each cusp has a preferred marking. In the case of a torus cusp, this is pair of essential simple curves meeting in one point; equivalently, a basis of the first homology of the boundary torus. These curves are called the meridian and the longitude.

This method changes these markings in various ways. In many cases, if the flag return_matrices is True then it returns a list of change-of-basis matrices is returned, one per cusp, which will restore the original markings if passed as peripheral_data.

simplify(passes_at_fours=6)

Try to simplify the triangulation by doing Pachner moves.

>>> M = Triangulation('12n123')
>>> M.simplify()

It does four kinds of moves that reduce the number of tetrahedra:

  • 3 -> 2 and 2 -> 0 Pacher moves, which eliminate one or two tetrahedra respectively.

  • On suitable valence-1 edges, does a 2 -> 3 and then 2 -> 0 move, which removes a tetrahedron and creates a new valence-1 edge.

  • When a 2-simplex has two edges of valence-4 giving rise to the suspension of a pentagon, replace these 6 tetrahedra with a single edge of valence 5.

It also does random 4 -> 4 moves in hopes of setting up a simplfication. The argument passes_at_fours is the number of times it goes through the valence-4 edges without progress before giving up.

slice_obstruction_HKL(primes_spec, verbose=False, check_in_S3=True)

For the exterior of a knot in S^3, searches for a topological slicing obstruction from:

Herald, Kirk, Livingston, Math Zeit., 2010 https://dx.doi.org/10.1007/s00209-009-0548-1 https://arxiv.org/abs/0804.1355

The test looks at the cyclic branched covers of the knot of prime order p and the F_q homology thereof where q is an odd prime. The range of such (p, q) pairs searched is given by primes_spec as a list of (p_max, [q_min, q_max]). It returns the pair (p, q) of the first nonzero obstruction found (in which case K is not slice), and otherwise returns None:

sage: M = Manifold('K12n813')
sage: spec = [(10, [0, 20]), (20, [0, 10])]
sage: M.slice_obstruction_HKL(spec, verbose=True)
   Looking at (2, 3) ...
   Looking at (3, 2) ...
   Looking at (3, 7) ...
(3, 7)

You can also specify the p to examine by a range [p_min, p_max] or the q by just q_max:

sage: spec = [([3, 10], 10)]
sage: M.slice_obstruction_HKL(spec, verbose=True)
   Looking at (3, 2) ...
   Looking at (3, 7) ...
(3, 7)

If primes_spec is just a pair (p, q) then only that obstruction is checked:

sage: M.slice_obstruction_HKL((2, 3))
sage: M.slice_obstruction_HKL((3, 7))
(3, 7)

Technical note: As implemented, can only get an obstruction when the decomposition of H_1(cover; F_q) into irreducible Z/pZ-modules has no repeat factors. The method of [HKL] can be used more broadly, but other cases requires computing many more twisted Alexander polynomials.

triangulation_isosig(decorated=True, ignore_cusp_ordering=False, ignore_curve_orientations=False, ignore_orientation=True)

Returns a compact text representation of the triangulation, called a “decorated isomorphism signature”

>>> T = Triangulation('m004')
>>> T.triangulation_isosig()
'cPcbbbiht_BaCB'

You can use this string to recreate an isomorphic triangulation later

>>> A = Triangulation('y233')
>>> A.triangulation_isosig()
'hLMzMkbcdefggghhhqxqhx_BaaB'
>>> B = Triangulation('hLMzMkbcdefggghhhqxqhx_BaaB')
>>> A == B
True

By default, the returned string encodes the peripheral curves (and slopes of Dehn-fillings if any are present), but you can request only the “isomorphism signature” which can be given to Regina.

>>> E = Triangulation('K3_1')   # the (-2, 3, 7) exterior
>>> isosig = E.triangulation_isosig(decorated = False); isosig
'dLQacccjsnk'
>>> F = Triangulation(isosig)
>>> E.isomorphisms_to(F)[1]
0 -> 0
[1 18]
[0  1]
Extends to link
>>> E.triangulation_isosig()
'dLQacccjsnk_BaRsB'
>>> F.triangulation_isosig()
'dLQacccjsnk_BaaB'
>>> G = Triangulation('dLQacccjsnk_BaRsB')
>>> E.isomorphisms_to(G)[0]
0 -> 0
[1 0]
[0 1]
Extends to link

If you do not care about the indexing of the cusps when using a decorated signature, use ignore_cusp_ordering

>>> M = Manifold("L14n64110(1,2)(2,3)(-2,1)(3,4)(0,0)")
>>> isosig = M.triangulation_isosig(decorated = True, ignore_cusp_ordering = True)
>>> isosig
'xLLvLvMLPMPLAMQQcceflnjmmmospsrttvvvtswwwiieiifdeauinasltltahmbjn_bacBbaaBBaBbBbbaabba(2,3)(-2,1)(1,2)(3,4)(0,0)'
>>> N = Manifold(isosig).filled_triangulation()
>>> N.is_isometric_to(M.filled_triangulation())
True

If you do not care about the orientations of the peripheral curves, use ignore_curve_orientations

>>> M = Manifold("L6a1")
>>> M.triangulation_isosig()
'gLLAQcdeefffdopuado_BabbBaab'
>>> isosig = M.triangulation_isosig(decorated = True, ignore_curve_orientations = True)
>>> isosig
'gLLAQcdeefffdopuado_babbbaab'
>>> N = Manifold(isosig)
>>> M.isomorphisms_to(N)
[0 -> 0  1 -> 1
[-1 0]  [-1 0]
[ 0 1]  [ 0 1]
Extends to link, 0 -> 0  1 -> 1
[1  0]  [1  0]
[0 -1]  [0 -1]
Extends to link]

By default, the isomorphism signature does not capture the orientation of an orientable triangulation. If you specify ignore_orientation = False, the isomorphism signature for an oriented triangulation and its mirror image will be different if the triangulation is cheiral.

>>> M = Manifold("m006")
>>> M.triangulation_isosig(decorated = False, ignore_orientation = False)
'dLQacccjnjs'
>>> M.reverse_orientation()
>>> M.triangulation_isosig(decorated = False, ignore_orientation = False)
'dLQacccnsnk'

Note that a decorated triangulation isosig with the default values ignore_orientation = True but ignore_curve_orientations = False still captures the orientations of the triangulation through the peripheral curves.

>>> M = Manifold("m006")
>>> M.triangulation_isosig()
'dLQacccjnjs_aBbB'
>>> M.reverse_orientation()
>>> M.triangulation_isosig()
'dLQacccjnjs_aBBb'

The code has been copied from Regina where the corresponding method is called “isoSig”.

Unlike dehydrations for 3-manifold triangulations, an isomorphism signature uniquely determines a triangulation up to combinatorial isomorphism. That is, two triangulations of 3-dimensional manifolds are combinatorially isomorphic if and only if their isomorphism signatures are the same string. For full details, see Simplification paths in the Pachner graphs of closed orientable 3-manifold triangulations, Burton, 2011.

For details about how the peripheral decorations work, see the SnapPy source code.

with_hyperbolic_structure()

Add a (possibly degenerate) hyperbolic structure, turning the Triangulation into a Manifold.

>>> M = Triangulation('m004')
>>> N = M.with_hyperbolic_structure()
>>> N.volume() 
2.02988321