Manifold: the main class

class snappy.Manifold(spec=None)

A Manifold is a Triangulation together with a geometric structure. That is, a Manifold is an ideal triangulation of the interior of a compact 3-manifold with torus boundary components, where each tetrahedron has been assigned the geometry of an ideal tetrahedron in hyperbolic 3-space. A Dehn-filling can be specified for each boundary component, allowing the description of closed 3-manifolds, some orbifolds and cone 3-manifolds. Here’s a quick example:

>>> M = Manifold('9_42')
>>> M.volume()  
4.05686022
>>> M.cusp_info('shape') 
[-4.278936315 + 1.95728679*I]

This is an example for running SnapPy inside Sage:

sage: import snappy
sage: M = snappy.Manifold('m125(1,2)(4,5)')
sage: M.is_orientable()
True

An alternative way of running SnapPy inside Sage:

sage: from snappy import *
sage: M = Manifold('m123')
sage: M.num_cusps()
1

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

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

  • Manifold(‘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).

  • Manifold() : 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’.

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

    • Callahan-Dean-Weeks-Champanerkar-Kofman-Patterson knots: e.g. ‘K6_21’.

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

  • 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. See Manifold.save() for details.

  • 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.

browse()

Opens browser window with a graphical interface, which allows to explore the manifold and interact with it. This includes: invariants, Dirichlet domain, cusp neighborhoods, inside view, symmetry, Dehn filling, drilling, etc.

>>> M = Manifold('m125')
>>> M.browse()  

This does not work when using SnapPy in a Docker container.

canonical_retriangulation(verified=False, interval_bits_precs=[53, 212], exact_bits_prec_and_degrees=[(212, 10), (1000, 20), (2000, 20)], verbose=False)

The canonical retriangulation which is closely related to the canonical cell decomposition and described in more detail here:

>>> M = Manifold("m412")
>>> K = M.canonical_retriangulation()
>>> len(K.isomorphisms_to(K)) # Unverified size of the isometry group.
8

When used inside Sage and verified = True is passed as argument, the verify module will certify the result to be correct:

sage: M = Manifold("m412")
sage: K = M.canonical_retriangulation(verified = True)
sage: len(K.isomorphisms_to(K)) # Verified size of the isometry group.
8

See verify.verified_canonical_retriangulation() for the additional options.

canonize()

Change the triangulation to an arbitrary retriangulation of the canonical cell decomposition.

>>> M = Manifold('m007')
>>> M.num_tetrahedra()
3
>>> M.canonize()
>>> M.num_tetrahedra()
4

Note: due to rounding error, it is possible that this is not actually the canonical triangulation.

chern_simons(accuracy=False)

Returns the Chern-Simons invariant of the manifold, if it is known.

>>> M = Manifold('m015')
>>> M.chern_simons() 
-0.15320413

The return value has an extra attribute, accuracy, which is the number of digits of accuracy as estimated by SnapPea.

>>> cs, accuracy = M.chern_simons(accuracy = True)
>>> accuracy in (8, 9, 56) # Low and High precision
True

By default, when the manifold has at least one cusp, Zickert’s algorithm is used; when the manifold is closed we use SnapPea’s original algorithm, which is based on Meyerhoff-Hodgson-Neumann.

Note: When computing the Chern-Simons invariant of a closed manifold, one must sometimes compute it first for the unfilled manifold so as to initialize SnapPea’s internals. For instance,

>>> M = Manifold('5_2')
>>> M.chern_simons() 
-0.15320413
>>> M.dehn_fill( (1,2) )
>>> M.chern_simons() 
0.07731787

works, but will fail with ‘Chern-Simons invariant not currently known’ if the first call to chern_simons is not made.

complex_volume(verified_modulo_2_torsion=False, bits_prec=None)

Returns the complex volume, i.e. volume + i 2 pi^2 (chern simons)

>>> M = Manifold('5_2')
>>> M.complex_volume() 
2.82812209 - 3.02412838*I
>>> c = M.chern_simons()
>>> M.dehn_fill((1,2))
>>> M.complex_volume() 
2.22671790 + 1.52619361*I
>>> M = Manifold("3_1")
>>> cvol = M.complex_volume()
>>> cvol.real() 
0
>>> cvol.imag() 
-1.64493407

If no cusp is filled or there is only one cusped (filled or unfilled), the complex volume can be verified up to multiples of i pi^2 / 2 by passing verified_modulo_2_torsion = True when inside SageMath (and higher precision can be requested with bits_prec):

sage: M = Manifold("m015")
sage: M.complex_volume(verified_modulo_2_torsion=True, bits_prec = 90) # doctest: +NUMERIC24
2.828122088330783162764? + 1.910673824035377649698?*I
sage: M = Manifold("m015(3,4)")
sage: M.complex_volume(verified_modulo_2_torsion=True) # doctest: +NUMERIC6
2.625051576? - 0.537092383?*I
copy()

Returns a copy of the manifold

>>> M = Manifold('m015')
>>> N = M.copy()
cover(permutation_rep)

M.cover(permutation_rep)

Returns a Manifold 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 = Manifold('m004')
>>> N0 = M.cover([[1, 3, 0, 4, 2], [0, 2, 1, 4, 3]])
>>> abs(N0.volume()/M.volume() - 5) < 0.0000000001
True

If 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. Some examples:

sage: M = Manifold('m004')

The basic method:

sage: N0 = M.cover([[1, 3, 0, 4, 2], [0, 2, 1, 4, 3]])

From a Gap subgroup:

sage: G = gap(M.fundamental_group())
sage: H = G.LowIndexSubgroupsFpGroup(5)[9]
sage: N1 = M.cover(H)
sage: N0 == N1
True

Or a homomorphism to a permutation group:

sage: f = G.GQuotients(PSL(2,7))[1]
sage: N2 = M.cover(f)
sage: N2.volume()/M.volume() # doctest: +NUMERIC9
8.00000000

Or maybe we want larger cover coming from the kernel of this:

sage: N3 = M.cover(f.Kernel())
sage: N3.volume()/M.volume() # doctest: +NUMERIC9
168.00000000

Check the homology against what Gap computes directly:

sage: N3.homology().betti_number()
32
sage: len([ x for x in f.Kernel().AbelianInvariants().sage() if x == 0])
32

We can do the same for Magma:

sage: G = magma(M.fundamental_group())             #doctest: +SKIP
sage: Q, f = G.pQuotient(5, 1, nvals = 2)          #doctest: +SKIP
sage: M.cover(f.Kernel()).volume()                 #doctest: +SKIP
10.14941606
sage: h = G.SimpleQuotients(1, 11, 2, 10000)[1,1]  #doctest: +SKIP
sage: N4 = M.cover(h)                              #doctest: +SKIP
sage: N2 == N4                                     #doctest: +SKIP
True
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 Manifolds 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 = Manifold('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 = Manifold('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_area_matrix(method='trigDependentTryCanonize', verified=False, bits_prec=None)

This function returns a matrix that can be used to check whether cusp neighborhoods of areas a0, …, am-1 are disjoint: the cusp neighborhoods about cusp i and j are disjoint (respectively, the cusp neighborhood embeds if i and j are equal) if ai * aj is less than or equal to the entry (i,j) of the cusp area matrix. Note that the “if” becomes “if and only if” if we pick the “maximal cusp area matrix”.

This function can operate in different ways (determined by method). By default (method='trigDependentTryCanonize'), it returns a result which can be suboptimal and non-deterministic but is quicker to compute and sufficies for many applications:

>>> M = Manifold("s776")
>>> M.cusp_area_matrix() 
[28.0000000000000 7.00000000000000 6.99999999999999]
[7.00000000000000 21.4375000000000 7.00000000000000]
[6.99999999999999 7.00000000000000 21.4375000000000]

If method='maximal' is specified, the result is the “maximal cusp area matrix”, thus it is optimal and an invariant of the manifold with labeled cusps. Note that the “maximal cusp area matrix” is only available as verified computation and thus requires passing verified = True:

sage: M.cusp_area_matrix(method = 'maximal', verified=True) # doctest: +NUMERIC6
[28.0000000000?  7.0000000000?  7.0000000000?]
[ 7.0000000000?  28.000000000? 7.00000000000?]
[ 7.0000000000? 7.00000000000?   28.00000000?]

If verified = True is specified and method is not maximal, the entries are all guaranteed to be less than the corresponding ones in the maximal cusp area matrix (more precisely, the lower end point of the interval is guaranteed to be less than the true value of the corresponding maximal cusp area matrix entry):

sage: M.cusp_area_matrix(verified=True, bits_prec=70) # doctest: +NUMERIC15
[ 28.000000000000000?  7.0000000000000000?  7.0000000000000000?]
[ 7.0000000000000000? 21.4375000000000000?  7.0000000000000000?]
[ 7.0000000000000000?  7.0000000000000000? 21.4375000000000000?]

For expert users:

Besides the two values above, method can be trigDependent: this result is also fast to compute by making the assumption that cusp neighborhoods are not only disjoint but also in “standard form” with respect to the triangulation (i.e., when lifting of a cusp neighborhood to a horoball in the universal cover, it intersects a geodesic tetrahedron in three but not four faces). trigDependentTryCanonize is similar to trigDependent but tries to “proto-canonize” (a copy of) the triangulation first since this often produces a matrix that is closer to the maximal cusp area matrix, for example:

>>> M = Manifold("o9_35953")
>>> M.cusp_area_matrix(method = 'trigDependent') 
[72.9848715318467 12.7560424258060]
[12.7560424258060 6.65567118002656]
>>> M.cusp_area_matrix(method = 'trigDependentTryCanonize') 
[72.9848715318466 12.7560424258060]
[12.7560424258060 62.1043047674605]

Compare to maximal area matrix:

sage: M.cusp_area_matrix(method = 'maximal', verified = True, bits_prec = 100) # doctest: +NUMERIC15
[       72.984871531846664? 12.7560424258059765562778?]
[12.7560424258059765562778?     62.104304767460978078?]
cusp_areas(policy='unbiased', method='trigDependentTryCanonize', verified=False, bits_prec=None, first_cusps=[])

Picks areas for the cusps such that the corresponding cusp neighborhoods are disjoint. By default, the policy is unbiased which means that the cusp neighborhoods are blown up simultaneously with a cusp neighborhood stopping to grow when it touches another cusp neighborhood or itself:

>>> M = Manifold("s776")
>>> M.cusp_areas() 
[2.64575131106459, 2.64575131106459, 2.64575131106459]

Alternatively, policy='greedy' means that the first cusp neighborhood is blown up until it touches itself, then the second cusp neighborhood is blown up until it touches itself or the first cusp neighborhood, …:

>>> M.cusp_areas(policy='greedy') 
[5.29150262212918, 1.32287565553230, 1.32287565553229]

To specify cusps to be blown up first, and in which order to blow them up, set first_cusps to the appropriate list of cusps.

>>> M = Manifold('o9_44210')
>>> M.cusp_areas(policy='greedy') 
[7.053940530873898, 3.2712450270, 2.7091590087]
>>> M.cusp_areas(policy='greedy', first_cusps=[]) 
[7.053940530873898, 3.2712450270, 2.7091590087]
>>> M.cusp_areas(policy='greedy', first_cusps=[0,]) 
[7.053940530873898, 3.2712450270, 2.7091590087]
>>> M.cusp_areas(policy='greedy', first_cusps=[0,1]) 
[7.053940530873898, 3.2712450270, 2.7091590087]
>>> M.cusp_areas(policy='greedy', first_cusps=[0,1,2]) 
[7.053940530873898, 3.2712450270, 2.7091590087]
>>> M.cusp_areas(policy='greedy', first_cusps=[0,2,1]) 
[7.053940530873898, 2.3513135103, 3.7690945490]
>>> M.cusp_areas(policy='greedy', first_cusps=[1,]) 
[4.0302253322, 5.725527974287718, 1.5478612583]

cusp_areas is implemented using Manifold.cusp_area_matrix() and the same arguments (method, verified, bits_prec) are accepted. For example, verified computations are supported:

sage: M=Manifold("v2854")
sage: M.cusp_areas(verified=True) # doctest: +NUMERIC9
[3.6005032476?, 3.6005032476?]

If method='maximal', policy='unbiased' and verified=True, the result is an invariant of the manifold with labeled cusps and the corresponding cusp neighborhoods are maximal in that every cusp neighborhood is touching some (not necessarily distinct) cusp neighborhood.

Area of the cusp neighborhood touching itself for a one-cusped manifold:

sage: M=Manifold("v1959")
sage: M.cusp_areas(method='maximal', verified=True, bits_prec=100) # doctest: +NUMERIC15
[7.15679216175810579?]
cusp_info(data_spec=None, verified=False, bits_prec=None)

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

>>> M = Manifold('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)

To get more detailed information about the cusp, we do

>>> c = M.cusp_info(0)
>>> c.shape 
0.11044502 + 0.94677098*I
>>> c.modulus
-0.12155872 + 1.04204128*I
>>> sorted(c.keys())
['filling', 'holonomies', 'holonomy_accuracy', 'index', 'is_complete', 'modulus', 'shape', 'shape_accuracy', 'topology']

Here ‘shape’ is the shape of the cusp, i.e. (longitude/meridian) and ‘modulus’ is its shape in the geometrically preferred basis, i.e. ( (second shortest translation)/(shortest translation)). For cusps that are filled, one instead cares about the holonomies:

>>> M.cusp_info(-1)['holonomies'] 
(-0.59883089 + 1.09812548*I, 0.89824633 + 1.49440443*I)

The complex numbers returned for the shape and for the two holonomies have an extra attribute, accuracy, which is SnapPea’s estimate of their accuracy.

You can also get information about multiple cusps at once:

>>> M.cusp_info() 
[Cusp 0 : complete torus cusp of shape 0.11044502 + 0.94677098*I,
 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]

The cusp shapes can be verified:

sage: M = Manifold('m292')
sage: M.cusp_info('shape', verified = True, bits_prec = 60) # doctest: +NUMERIC12
[-0.1766049820997? + 1.2028208192855?*I,
 -0.1766049820997? + 1.2028208192855?*I]
cusp_neighborhood()

Returns information about the cusp neighborhoods of the manifold, in the form of data about the corresponding horoball diagrams in hyperbolic 3-space.

>>> M = Manifold('s000')
>>> CN = M.cusp_neighborhood()
>>> CN.volume() 
0.32475953
>>> len(CN.horoballs(0.01))
178
>>> CN.view()  # Opens picture of the horoballs  
cusp_translations(policy='unbiased', method='trigDependentTryCanonize', verified=False, bits_prec=None, first_cusps=[])

Picks disjoint cusp neighborhoods and returns the respective (complex) Euclidean translations of the meridian and longitude for each cusp. The method takes the same arguments as Manifold.cusp_areas() and uses that method to pick the cusp neighborhood. The result is a list of pairs, the second entry corresponding to a longitude is always real:

>>> M = Manifold("s776")
>>> M.cusp_translations() 
[(0.500000000000000 + 1.32287565553230*I, 2.00000000000000), (0.500000000000000 + 1.32287565553230*I, 2.00000000000000), (0.499999999999999 + 1.32287565553230*I, 2.00000000000000)]

Arguments such as policy='greedy' are interpreted the same way as for Manifold.cusp_areas():

>>> M.cusp_translations(policy = 'greedy', first_cusps = [], bits_prec = 100) 
[(0.70710678118654752440084436210 + 1.8708286933869706927918743662*I, 2.8284271247461900976033774484), (0.35355339059327376220042218105 + 0.93541434669348534639593718308*I, 1.4142135623730950488016887242), (0.35355339059327376220042218105 + 0.93541434669348534639593718308*I, 1.4142135623730950488016887242)]

and can return verified intervals:

sage: M.cusp_translations(method = 'maximal', verified = True) # doctest: +NUMERIC9
[(0.50000000000? + 1.32287565553?*I, 2.00000000000?), (0.500000000000? + 1.32287565554?*I, 2.00000000000?), (0.500000000000? + 1.32287565554?*I, 2.00000000000?)]

that are guaranteed to contain the true translations of cusp neighborhoods verified to be disjoint (the element corresponding to a longitude is always in a RealIntervalField).

Remark: The default method = 'trigDependentTryCanonize' is (potentially) non-deterministic and thus the result of

[ M.cusp_translations()[i] for i in range(M.num_cusps()) ]

might not correspond to disjoint cusp neighborhoods.

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 = Manifold('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 = Manifold('m004')
    >>> N.dehn_fill( (-3,4) )
    >>> N
    m004(-3,4)
    

Does not return a new Manifold.

dirichlet_domain(vertex_epsilon=default_vertex_epsilon, displacement=(0.0, 0.0, 0.0), centroid_at_origin=True, maximize_injectivity_radius=True, include_words=False)

Returns a DirichletDomain object representing a Dirichlet domain of the hyperbolic manifold, typically centered at a point which is a local maximum of injectivity radius. It will have ideal vertices if the manifold is not closed.

>>> M = Manifold('m015')
>>> D = M.dirichlet_domain()
>>> D
32 finite vertices, 2 ideal vertices; 54 edges; 22 faces
>>> D.view()   #Shows 3d-graphical view.  

The group elements for the face-pairings of the Dirichlet domain can be given as words in the original generators of the (unsimplified) fundamental group by setting include_words = True:

>>> sorted(M.dirichlet_domain(include_words = True).pairing_words()) 
['A', ...]

Other options can be provided to customize the computation; the default choices are shown below:

>>> M.dirichlet_domain(vertex_epsilon=10.0**-8,
...    displacement = [0.0, 0.0, 0.0],
...    centroid_at_origin=True, maximize_injectivity_radius=True)
32 finite vertices, 2 ideal vertices; 54 edges; 22 faces

Here’s one with different combinatorics:

>>> E = M.dirichlet_domain(displacement=[0.5, 0.3, -0.2],
...                        maximize_injectivity_radius=False)
>>> E
44 finite vertices, 1 ideal vertices; 69 edges; 26 faces
drill(which_curve, max_segments=6)

Drills out the specified dual curve from among all dual curves with at most max_segments, which defaults to 6. The method dual_curve allows one to see the properties of curves before choosing which one to drill out.

>>> M = Manifold('v3000')
>>> N = M.drill(0, max_segments=3)
>>> (M.num_cusps(), N.num_cusps())
(1, 2)
drill_word(word: str, verified: bool = False, bits_prec=None, verbose: bool = False)

Drills the geodesic corresponding to the given word in the unsimplified fundamental group.

>>> from snappy import Manifold
>>> M = Manifold("m004")
>>> M.length_spectrum(1.2, include_words = True, grouped = False)
mult  length                                  topology     parity word
1     1.08707014499574 -  1.72276844987009*I  circle       +      a
1     1.08707014499574 +  1.72276844987009*I  circle       +      bC
>>> N = M.drill_word('a')
>>> N.identify()
[m129(0,0)(0,0), 5^2_1(0,0)(0,0), L5a1(0,0)(0,0), ooct01_00001(0,0)(0,0)]

The last cusp of the new manifold corresponds to the drilled geodesic and the longitude and meridian for that cusp are chosen such that (1,0)-filling results in the original (undrilled) manifold. The orientation of the new longitude is chosen so that it is parallel to the closed geodesic. That is, the new longitude is homotopic to the closed geodesic when embedding the drilled manifold into the original manifold.

>>> N.dehn_fill((1,0),1)
>>> M.is_isometric_to(N)
True
>>> N.cusp_info(1)['core_length'] 
1.08707014499574 - 1.72276844987009*I

If the drilled geodesic coincides with a core curve of a filled cusp, the cusp is unfilled instead and the longitude and meridian changed so that the above again applies. The cusp order is also changed so that the unfilled cusp becomes the last cusp.

>>> M = Manifold("m004(2,3)")
>>> M.volume() 
1.73712388065
>>> M.cusp_info(0)['core_length'] 
0.178792491242577 - 2.11983007979743*I
>>> M.fundamental_group(simplify_presentation = False).complex_length('aBAbbABab') 
0.178792491242577 - 2.11983007979743*I
>>> N = M.drill_word('aBAbbABab')
>>> N
m004_drilled(0,0)
>>> N.num_cusps()
1
>>> N.dehn_fill((1,0))
>>> N.volume() 
1.73712388065

Even though the output of the drilling geodesic algorithm is a triangulation and thus combinatorial in nature, the intermediate computations to compute the intersections of the geodesic with the faces of the tetrahedra is numerical. Sometimes it is necessary to increase the precision with bits_prec to make this computation accurate or succeed. If verified = True is specified, intervals are used for all computations and the result is provably correct (only supported when used inside SageMath). That is, the algorithm will fail with an exception (most likely InsufficientPrecisionError) if insufficient precision is used. Example of verified computation:

sage: M = Manifold("m004(2,3)")
sage: M.drill_word('caa', verified = True, bits_prec = 100)
m004_drilled(2,3)(0,0)

An example where we drill the core geodesic:

>>> from snappy import Manifold
>>> M = Manifold("v2986(3,4)")
>>> N = M.drill_word('EdFgabcGEdFgaDcc')
>>> N.is_isometric_to(Manifold("v2986"), return_isometries = True) 
[0 -> 0
 [3 -1]
 [4 -1]
 Does not extend to link]
drill_words(words: Sequence[str], verified: bool = False, bits_prec=None, verbose: bool = False)

A generalization of M.drill_word taking a list of words to drill several geodesics simultaneously.

Here is an example where we drill the core curve corresponding to the third cusp and a geodesic that is not a core curve:

>>> from snappy import Manifold
>>> M=Manifold("t12047(0,0)(1,3)(1,4)(1,5)")
>>> [ info.get('core_length') for info in M.cusp_info() ] 
[None,
 0.510804267610103 + 1.92397456664239*I,
 0.317363079597924 + 1.48157893409218*I,
 0.223574975263386 + 1.26933288854145*I]
>>> G = M.fundamental_group(simplify_presentation = False)
>>> G.complex_length('c') 
0.317363079597924 + 1.48157893409218*I
>>> G.complex_length('fA') 
1.43914411734250 + 2.66246879992795*I
>>> N = M.drill_words(['c','fA'])
>>> N
t12047_drilled(0,0)(1,3)(1,5)(0,0)(0,0)

The last n cusps correspond to the n geodesics that were drilled, appearing in the same order the words for the geodesics were given. Note that in the above example, the drilled manifold has only five cusps even though the original manifold had four cusps and we drilled two geodesics. This is because one geodesic was a core curve. The corresponding cusp was unfilled (from (1,4)) and grouped with the other cusps coming from drilling.

We obtain the original (undrilled) manifold by (1,0)-filling the last n cusps.

>>> N.dehn_fill((1,0), 3)
>>> N.dehn_fill((1,0), 4)
>>> M.is_isometric_to(N)
True
>>> [ info.get('core_length') for info in N.cusp_info() ] 
[None,
 0.510804267610103 + 1.92397456664239*I,
 0.223574975263386 + 1.26933288854145*I,
 0.317363079597924 + 1.48157893409218*I,
 1.43914411734251 + 2.66246879992796*I]
dual_curves(max_segments=6)

Constructs a reasonable selection of simple closed curves in a manifold’s dual 1-skeleton. In particular, it returns those that appear to represent geodesics. The resulting curves can be drilled out.

>>> M = Manifold('m015')
>>> curves = M.dual_curves()
>>> curves 
[  0: orientation-preserving curve of length 0.56239915 - 2.81543089*I,
   1: orientation-preserving curve of length 1.12479830 + 0.65232354*I,
   2: orientation-preserving curve of length 1.26080402 + 1.97804689*I,
   3: orientation-preserving curve of length 1.58826933 + 1.67347167*I,
   4: orientation-preserving curve of length 1.68719745 + 2.81543089*I]

Each curve is returned as an info object with these keys

>>> sorted(curves[0].keys())
['complete_length', 'filled_length', 'index', 'max_segments', 'parity']

We can drill out any of these curves to get a new manifold with one more cusp.

>>> N = M.drill(curves[0])
>>> (M.num_cusps(), N.num_cusps())
(1, 2)

By default, this function only finds curves of length 6; this can be changed by specifying the optional argument max_segments

>>> M.dual_curves(max_segments=2) 
[  0: orientation-preserving curve of length 0.56239915 - 2.81543089*I]
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.

Filling all the cusps results in a Triangulation rather than a Manifold, since SnapPea can’t deal with hyperbolic structures when there are no cusps.

Examples:

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

Filling cusps 0 and 2 :

>>> M = Manifold('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)

Return a HolonomyGroup representing the fundamental group of the manifold, together with its holonomy representation. If integer Dehn surgery parameters have been set, then the corresponding peripheral elements are killed.

>>> M = Manifold('m004')
>>> G = M.fundamental_group()
>>> G
Generators:
   a,b
Relators:
   aaabABBAb
>>> G.peripheral_curves()
[('ab', 'aBAbABab')]
>>> G.SL2C('baaBA') 
[ 2.50000000 - 2.59807621*I -6.06217783 - 0.50000000*I]
[ 0.86602540 - 2.50000000*I -4.00000000 + 1.73205081*I]

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

  • simplify_presentation

  • fillings_may_affect_generators

  • minimize_number_of_generators

>>> 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
high_precision()

Return a high precision version of this manifold.

>>> M = Manifold('m004')
>>> type(M.high_precision())
<class 'snappy.ManifoldHP'>
holonomy_matrix_entries(fundamental_group_args=[], match_kernel=True)

The entries of the matrices of the holonomy as list of ApproximateAlgebraicNumbers (four consecutive numbers per matrix). The numbers are guaranteed to lie in the trace field only if match_kernel = False:

sage: M = Manifold("m004")
sage: mat_entries = M.holonomy_matrix_entries(match_kernel = False) # doctest: +NORMALIZE_WHITESPACE +NUMERIC9
sage: mat_entries
<SetOfAAN: [0.5 + 0.8660254037844386*I, 0.5 - 0.8660254037844386*I, 0.5 + 0.8660254037844386*I, 1.0 - 1.7320508075688772*I, 1.0 - 3.4641016151377544*I, -2.0 + 1.7320508075688772*I, -1.0 - 1.7320508075688772*I, 1.7320508075688772*I]>
sage: K = mat_entries.find_field(100, 10, optimize = True)[0]
sage: K.polynomial()
x^2 - x + 1
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
hyperbolic_SLN_torsion(N, bits_prec=100)

Compute the torsion polynomial of the holonomy representation lifted to SL(2, C) and then followed by the irreducible representation from SL(2, C) -> SL(N, C):

sage: M = Manifold('m016')
sage: [M.hyperbolic_SLN_torsion(N).degree() for N in [2, 3, 4]]
[18, 27, 36]
hyperbolic_adjoint_torsion(bits_prec=100)

Computes the torsion polynomial of the adjoint representation a la Dubois-Yamaguichi. This is not a sign-refined computation so the result is only defined up to sign, not to mention a power of the variable ‘a’:

sage: M = Manifold('K11n42')
sage: tau = M.hyperbolic_adjoint_torsion()
sage: tau.parent()
Univariate Polynomial Ring in a over Complex Field with 100 bits of precision
sage: tau.degree()
7
hyperbolic_torsion(bits_prec=100, all_lifts=False, wada_conventions=False, phi=None)

Computes the hyperbolic torsion polynomial as defined in [DFJ]:

sage: M = Manifold('K11n42')
sage: M.alexander_polynomial()
1
sage: tau = M.hyperbolic_torsion(bits_prec=200)
sage: tau.degree()
6
identify(extends_to_link=False)

Looks for the manifold in all of the SnapPy databases. For hyperbolic manifolds this is done by searching for isometries:

>>> M = Manifold('m125')
>>> M.identify()
[m125(0,0)(0,0), L13n5885(0,0)(0,0), ooct01_00000(0,0)(0,0)]

By default, there is no restriction on the isometries. One can require that the isometry take meridians to meridians. This might return fewer results:

>>> M.identify(extends_to_link=True)
[m125(0,0)(0,0), ooct01_00000(0,0)(0,0)]

For closed manifolds, extends_to_link doesn’t make sense because of how the kernel code works:

>>> C = Manifold("m015(1,2)")
>>> C.identify()
[m006(-5,2)]
>>> C.identify(True)
[]
init_hyperbolic_structure(force_recompute=False)
inside_view(cohomology_class=None, geodesics=[])

Show raytraced inside view of hyperbolic manifold. See images and demo video.

>>> M = Manifold("m004")
>>> M.inside_view() 

Or show the cohomology fractal:

>>> M = Manifold("m004")
>>> M.inside_view(cohomology_class = 0) 

The cohomology class in H^2(M, bd M; R) producing the cohomology fractal can be specified as a cocycle or using an automatically computed basis (of, say, length n). Thus, cohomology_class can be one of the following.

  • An integer i between 0 and n - 1 to pick the i-th basis vector.

  • An array of length n specifying the cohomology class as linear combination of basis vectors.

  • A weight for each face of each tetrahedron.

invariant_trace_field_gens(fundamental_group_args=[])

The generators of the trace field as ApproximateAlgebraicNumbers. Can be used to compute the tetrahedra field, where the first two parameters are bits of precision and maximum degree of the field:

sage: M = Manifold('m007(3,1)')
sage: K = M.invariant_trace_field_gens().find_field(100, 10, optimize=True)[0]
sage: L = M.trace_field_gens().find_field(100, 10, optimize=True)[0]
sage: K.polynomial(), L.polynomial()
(x^2 - x + 1, x^4 - 2*x^3 + x^2 + 6*x + 3)
is_isometric_to(Manifold other, return_isometries=False)

Returns True if M and N are isometric, False if they not. A RuntimeError is raised in cases where the SnapPea kernel fails to determine either answer. (This is fairly common for closed manifolds.)

>>> M = Manifold('m004')
>>> N = Manifold('4_1')
>>> K = Manifold('5_2')
>>> M.is_isometric_to(N)
True
>>> N.is_isometric_to(K)
False

We can also get a complete list of isometries between the two manifolds:

>>> M = Manifold('5^2_1')  # The Whitehead link
>>> N = Manifold('m129')
>>> isoms = M.is_isometric_to(N, return_isometries = True)
>>> isoms[6]  # Includes action on cusps
0 -> 1  1 -> 0
[1  2]  [-1 -2]
[0 -1]  [ 0  1]
Extends 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.

Note: The answer True is rigorous, but the answer False may not be as there could be numerical errors resulting in finding an incorrect canonical triangulation.

is_orientable()

Return whether the underlying 3-manifold is orientable.

>>> M = Triangulation('x124')
>>> M.is_orientable()
False
is_two_bridge()

If the manifold is the complement of a two-bridge knot or link in S^3, then this method returns (p,q) where p/q is the fraction describing the link. Otherwise, returns False.

>>> M = Manifold('m004')
>>> M.is_two_bridge()
(2, 5)
>>> M = Manifold('m016')
>>> M.is_two_bridge()
False

Note: An answer of ‘True’ is rigorous, but not the answer ‘False’, as there could be numerical errors resulting in finding an incorrect canonical triangulation.

isometry_signature(of_link=False, verified=False, interval_bits_precs=[53, 212], exact_bits_prec_and_degrees=[(212, 10), (1000, 20), (2000, 20)], verbose=False)

The isomorphism signature of the canonical retriangulation. This is a complete invariant of the isometry type of a hyperbolic 3-manifold and described in more detail here:

>>> M = Manifold("m125")
>>> M.isometry_signature() # Unverified isometry signature
'gLLPQccdefffqffqqof'

When used inside Sage and verified = True is passed as argument, the verify module will certify the result to be correct:

sage: M = Manifold("m125")
sage: M.isometry_signature(verified = True) # Verified isometry signature
'gLLPQccdefffqffqqof'

When of_link = True is specified, the peripheral curves are included in such a way that the result is a complete invariant of a link. In particular, isometry_signature(of_link=True) is invariant under changing the ordering or orientations of the components or flipping all crossings of a link simultaneously (it passes ignore_cusp_order = True, ignore_curve_orientations = True to Manifold.triangulation_isosig()):

>>> Manifold("5^2_1").isometry_signature(of_link = True)
'eLPkbdcddhgggb_baCbbaCb'
>>> Manifold("7^2_8").isometry_signature(of_link = True)
'eLPkbdcddhgggb_bBcBbaCb'

See verify.verified_canonical_retriangulation() for the additional options.

Note that interval methods cannot verify a canonical retriangulation with non-tetrahedral cells such as in the cas of m412, so the following call returns None:

sage: M = Manifold("m412")
sage: M.isometry_signature(verified = True, exact_bits_prec_and_degrees = None)
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.

length_spectrum(cutoff=1.0, full_rigor=True, grouped=True, include_words=False)

Returns a list of geodesics (with multiplicities) of length up to the specified cutoff value. (The default cutoff is 1.0.)

Here’s a quick example:

>>> L = Manifold("m016").length_spectrum(0.75, include_words=True)
>>> L 
mult  length                                  topology     parity word
1     0.58460368501799 +  2.49537045556047*I  circle       +      a
1     0.72978937305180 +  3.02669828218116*I  circle       +      Bc

Access just the length:

>>> L[0].length 
0.584603685017987 + 2.495370455560469*I

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.
polished_holonomy(bits_prec=100, fundamental_group_args=[], lift_to_SL2=True, ignore_solution_type=False, dec_prec=None, match_kernel=True)

Return the fundamental group of M equipt with a high-precision version of the holonomy representation:

sage: M = Manifold('m004')
sage: G = M.polished_holonomy()
sage: G('a').trace()
1.5000000000000000000000000000 - 0.86602540378443864676372317075*I
sage: G = M.polished_holonomy(bits_prec=1000)
sage: G('a').trace().parent()
Complex Field with 1000 bits of precision
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.

  • Make the shortest curves the meridians, and the second shortest curves the longitudes.

    >>> M = Manifold('5_2')
    >>> M.cusp_info('shape') 
    [-2.49024467 + 2.97944707*I]
    >>> cob = M.set_peripheral_curves('shortest', return_matrices=True)
    >>> M.cusp_info('shape') 
    [-0.49024467 + 2.97944707*I]
    >>> cob
    [[[1, 0], [-2, 1]]]
    >>> M.set_peripheral_curves(cob)
    >>> M.cusp_info('shape') 
    [-2.49024467 + 2.97944707*I]
    

    You can also make just the meridians as short as possible while fixing the longitudes via the option ‘shortest_meridians’, and conversely with ‘shortest_longitudes’.

  • If cusps are Dehn filled, make those curves meridians.

    >>> M = Manifold('m125(0,0)(2,5)')
    >>> M.set_peripheral_curves('fillings')
    >>> M
    m125(0,0)(1,0)
    
  • Change the basis of a particular cusp, say the first one:

    >>> M.set_peripheral_curves( [ (1,2), (1,3) ] , 0)
    

    Here (1,2) is the new meridian written in the old basis, and (1,3) the new longitude.

  • Change the basis of all the cusps at once

    >>> new_curves = [ [(1,-1), (1,0)],  [(3,2), (-2,-1)] ]
    >>> M.set_peripheral_curves(new_curves)
    >>> M
    m125(0,0)(-1,-2)
    
set_target_holonomy(target, which_cusp=0, recompute=True)

M.set_target_holonomy(target, which_cusp=0, recompute=True)

Computes a geometric structure in which the Dehn filling curve on the specified cusp has holonomy equal to the target value. The holonomies of Dehn filling curves on other cusps are left unchanged. If the ‘recompute’ flag is False, the Dehn filling equations are modified, but not solved.

set_tetrahedra_shapes(filled_shapes=None, complete_shapes=None, fillings=None)

Replaces the tetrahedron shapes with those in the given lists, and sets the Dehn filling coefficients as specified by the fillings argument. The shapes will get double precision values; polishing will be needed for high precision shapes.

short_slopes(length=6, policy='unbiased', method='trigDependentTryCanonize', verified=False, bits_prec=None, first_cusps=[])

Picks disjoint cusp neighborhoods (using Manifold.cusp_areas(), thus the same arguments can be used) and returns for each cusp the slopes that have length less or equal to given length (defaults to 6) when measured on the boundary of the cusp neighborhood:

>>> M = Manifold("otet20_00022")
>>> M.short_slopes()
[[(1, 0), (-1, 1), (0, 1)], [(1, 0)]]

When verified=True, the result is guaranteed to contain all slopes of length less or equal to given length (and could contain additional slopes if precision is not high enough):

sage: M.short_slopes(verified = True)
[[(1, 0), (-1, 1), (0, 1)], [(1, 0)]]

The ten exceptional slopes of the figure-eight knot:

>>> M = Manifold("4_1")
>>> M.short_slopes()
[[(1, 0), (-4, 1), (-3, 1), (-2, 1), (-1, 1), (0, 1), (1, 1), (2, 1), (3, 1), (4, 1)]]

Two more slopes appear when increasing length to 2 pi:

>>> M.short_slopes(length = 6.283185307179586)
[[(1, 0), (-5, 1), (-4, 1), (-3, 1), (-2, 1), (-1, 1), (0, 1), (1, 1), (2, 1), (3, 1), (4, 1), (5, 1)]]

When using verified computations, length is converted into the RealIntervalField of requested precision:

sage: from sage.all import pi
sage: M.short_slopes(length = 2 * pi, verified = True, bits_prec = 100)
[[(1, 0), (-5, 1), (-4, 1), (-3, 1), (-2, 1), (-1, 1), (0, 1), (1, 1), (2, 1), (3, 1), (4, 1), (5, 1)]]
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.

solution_type(enum=False)

Returns the type of the current solution to the gluing equations, basically a summary of how degenerate the solution is. If the flag enum=True is set, then an integer value is returned. The possible answers are:

  • 0: ‘not attempted’

  • 1: ‘all tetrahedra positively oriented’ aka ‘geometric_solution’ Should correspond to a genuine hyperbolic structure.

  • 2: ‘contains negatively oriented tetrahedra’ aka ‘nongeometric_solution’ Probably corresponds to a hyperbolic structure but some simplices have reversed orientations.

  • 3: ‘contains flat tetrahedra’ All tetrahedra have shape in R - {0, 1}.

  • 4: ‘contains degenerate tetrahedra’ Some shapes are close to {0,1, or infinity}.

  • 5: ‘unrecognized solution type’

  • 6: ‘no solution found’

>>> M = Manifold('m007')
>>> M.solution_type()
'all tetrahedra positively oriented'
>>> M.dehn_fill( (3,1) )
>>> M.solution_type()
'contains negatively oriented tetrahedra'
>>> M.dehn_fill( (3,-1) )
>>> M.solution_type()
'contains degenerate tetrahedra'
split(which_surface)

Split the manifold open along a surface of positive characteristic found by the method “splitting_surfaces”. Returns a list of the pieces, with any sphere boundary components filled in.

Here’s an example of a Whitehead double on the trefoil.

>>> M = Manifold('K14n26039')
>>> S = M.splitting_surfaces()[0]
>>> S
Orientable two-sided with euler = 0
>>> pieces = M.split(S); pieces
[K14n26039.a(0,0)(0,0), K14n26039.b(0,0)]
>>> pieces[0].volume() 
3.66386238
>>> pieces[1].fundamental_group().relators()
['aabbb']

You can also specify a surface by its index.

>>> M = Manifold('L10n111')
>>> max( P.volume() for P in M.split(0) ) 
5.33348957
splitting_surfaces()

Searches for connected closed normal surfaces of nonnegative Euler characteristic. If spheres or projective planes are found, then tori and Klein bottles aren’t reported. There is no guarantee that all such normal surfaces will be found nor that any given surface is incompressible. The search is confined to surfaces whose quads are in the tetrahedra that have degenerate shapes.

You can split the manifold open along one of these surfaces using the method “split”.

A connect sum of two trefoils:

>>> M1 = Manifold('DT: fafBCAEFD')
>>> len(M1.splitting_surfaces())
2

First satellite knot in the table.

>>> M2 = Manifold('K13n4587')
>>> M2.splitting_surfaces()
[Orientable two-sided with euler = 0]
symmetric_triangulation()

Returns a Dehn filling description of the manifold realizing the symmetry group.

>>> M = Manifold('m003(-3,1)')
>>> M.symmetry_group()
D6
>>> N = M.symmetric_triangulation()
>>> N
m003(1,0)(1,0)(1,0)
>>> N.dehn_fill( [(0,0), (0,0), (0,0)] )
>>> N.symmetry_group(of_link=True)
D6
symmetry_group(of_link=False)

Returns the symmetry group of the Manifold. If the flag “of_link” is set, then it only returns symmetries that preserves the meridians.

tetrahedra_field_gens()

The shapes of the tetrahedra as ApproximateAlgebraicNumbers. Can be used to compute the tetrahedra field, where the first two parameters are bits of precision and maximum degree of the field:

sage: M = Manifold('m015')
sage: tets = M.tetrahedra_field_gens()
sage: tets.find_field(100, 10, optimize=True)    # doctest: +NORMALIZE_WHITESPACE +NUMERIC9
(Number Field in z with defining polynomial x^3 - x - 1
 with z = -0.6623589786223730? - 0.5622795120623013?*I,
<ApproxAN: -0.662358978622 - 0.562279512062*I>, [-z, -z, -z])
tetrahedra_shapes(part=None, fixed_alignment=True, bits_prec=None, dec_prec=None, intervals=False)

Gives the shapes of the tetrahedra in the current solution to the gluing equations. Returns a list containing one info object for each tetrahedron. The keys are:

  • rect : the shape of the tetrahedron, as a point in the complex plane.

  • log : the log of the shape

  • accuracies: a list of the approximate accuracies of the shapes, in order (rect re, rect im, log re, log im)

If the optional variable ‘part’ is set to one of the above, then the function returns only that component of the data.

If the flag ‘fixed_alignment’ is set to False, then the edges used to report the shape parameters are chosen so as to normalize the triangle.

>>> M = Manifold('m015')
>>> M.tetrahedra_shapes(part='rect') 
[0.66235898 + 0.56227951*I, 0.66235898 + 0.56227951*I, 0.66235898 + 0.56227951*I]
>>> M.tetrahedra_shapes() 
[{'accuracies': (11, 11, 12, 11), 'log': -0.14059979 + 0.70385772*I, 'rect': 0.66235898 + 0.56227951*I},
 {'accuracies': (11, 11, 11, 11), 'log': -0.14059979 + 0.70385772*I, 'rect': 0.66235898 + 0.56227951*I},
 {'accuracies': (11, 11, 11, 11), 'log': -0.14059979 + 0.70385772*I, 'rect': 0.66235898 + 0.56227951*I}]
trace_field_gens(fundamental_group_args=[])

The generators of the trace field as ApproximateAlgebraicNumbers. Can be used to compute the tetrahedra field, where the first two parameters are bits of precision and maximum degree of the field:

sage: M = Manifold('m125')
sage: traces = M.trace_field_gens()
sage: traces.find_field(100, 10, optimize=True)    # doctest: +NORMALIZE_WHITESPACE
(Number Field in z with defining polynomial x^2 + 1
 with z = -1*I,
<ApproxAN: -1.0*I>, [z + 1, z, z + 1])
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.

use_field_conversion(type cls, func)

A class method for specifying a numerical conversion function. This method is deprecated: SnapPy will automatically use SageMath number types or its own SnapPy number type depending on whether SageMath is available or not.

SnapPy includes its own number type, snappy.Number, which can represent floating point real or complex numbers of varying precision. (In fact, Number is a wrapper for a pari number of type ‘t_INT’, ‘t_FRAC’, ‘t_REAL’ or ‘t_COMPLEX’, and the pari gen can be extracted as an attribute: x.gen .) Methods of SnapPy objects which return numerical values will first compute the value as a Number, and then optionally convert the Number to a different numerical type which can be specified by calling this class method.

By default SnapPy returns Numbers when loaded into python, and elements of a Sage RealField or ComplexField when loaded into Sage. These will be 64 bit numbers for ordinary Manifolds and 212 bit numbers for high precision manifolds.

The func argument should be a function which accepts a number and returns a numerical type of your choosing. Alternatively, the strings ‘sage’ or ‘snappy’ can be passed as arguments to select either of the two default behaviors.

EXAMPLE:

sage: M = Manifold('m004')
sage: parent(M.volume())
Real Field with 64 bits of precision
sage: Manifold.use_field_conversion('snappy')
sage: M = Manifold('m004')
sage: parent(M.volume())
SnapPy Numbers with 64 bits precision
sage: Manifold.use_field_conversion('sage')
sage: M = Manifold('m004')
sage: parent(M.volume())
Real Field with 64 bits of precision
verify_hyperbolicity(verbose=False, bits_prec=None, holonomy=False, fundamental_group_args=[], lift_to_SL=True)

Given an orientable SnapPy Manifold, verifies its hyperbolicity.

Similar to HIKMOT’s verify_hyperbolicity(), the result is either (True, listOfShapeIntervals) or (False, []) if verification failed. listOfShapesIntervals is a list of complex intervals (elements in sage’s ComplexIntervalField) certified to contain the true shapes for the hyperbolic manifold.

Higher precision intervals can be obtained by setting bits_prec:

sage: from snappy import Manifold
sage: M = Manifold("m019")
sage: M.verify_hyperbolicity() # doctest: +NUMERIC12
(True, [0.780552527850? + 0.914473662967?*I, 0.780552527850? + 0.91447366296773?*I, 0.4600211755737? + 0.6326241936052?*I])

sage: M = Manifold("t02333(3,4)")
sage: M.verify_hyperbolicity() # doctest: +NUMERIC9
(True, [2.152188153612? + 0.284940667895?*I, 1.92308491369? + 1.10360701507?*I, 0.014388591584? + 0.143084469681?*I, -2.5493670288? + 3.7453498408?*I, 0.142120333822? + 0.176540027036?*I, 0.504866865874? + 0.82829881681?*I, 0.50479249917? + 0.98036162786?*I, -0.589495705074? + 0.81267480427?*I])

One can instead get a holonomy representation associated to the verified hyperbolic structure. This representation takes values in 2x2 matrices with entries in the ComplexIntervalField:

sage: M = Manifold("m004(1,2)")
sage: success, rho = M.verify_hyperbolicity(holonomy=True)
sage: success
True
sage: trace = rho('aaB').trace(); trace # doctest: +NUMERIC9
-0.1118628555? + 3.8536121048?*I
sage: (trace - 2).contains_zero()
False
sage: (rho('aBAbaabAB').trace() - 2).contains_zero()
True

Here, there is provably a fixed holonomy representation rho0 from the fundamental group G of M to SL(2, C) so that for each element g of G the matrix rho0(g) is contained in rho(g). In particular, the above constitutes a proof that the word ‘aaB’ is non-trivial in G. In contrast, the final computation is consistent with ‘aBAbaabAB’ being trivial in G, but does not prove this.

A non-hyperbolic manifold (False indicates that the manifold might not be hyperbolic but does not certify non-hyperbolicity. Sometimes, hyperbolicity can only be verified after increasing the precision):

sage: M = Manifold("4_1(1,0)")
sage: M.verify_hyperbolicity()
(False, [])

Under the hood, the function will call the CertifiedShapesEngine to produce intervals certified to contain a solution to the rectangular gluing equations. It then calls check_logarithmic_gluing_equations_and_positively_oriented_tets to verify that the logarithmic gluing equations are fulfilled and that all tetrahedra are positively oriented.

volume(accuracy=False, verified=False, bits_prec=None)

Returns the volume of the current solution to the hyperbolic gluing equations; if the solution is sufficiently non-degenerate, this is the sum of the volumes of the hyperbolic pieces in the geometric decomposition of the manifold.

>>> M = Manifold('m004')
>>> M.volume() 
2.02988321
>>> M.solution_type()
'all tetrahedra positively oriented'

The return value has an extra attribute, accuracy, which is the number of digits of accuracy as estimated by SnapPea. When printing the volume, the result is rounded to 1 more than this number of digits.

>>> vol, accuracy = M.volume(accuracy = True)
>>> accuracy in (10, 63) # Low precision, High precision
True

Inside SageMath, verified computation of the volume of a hyperbolic manifold is also possible (this will verify first that the manifold is indeed hyperbolic):

sage: M.volume(verified=True, bits_prec=100)   #doctest: +NUMERIC24
2.029883212819307250042405109?
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
without_hyperbolic_structure()

Returns self as a Triangulation, forgetting the hyperbolic structure in the process.

>>> M = Manifold('9_42')
>>> T = M.without_hyperbolic_structure()
>>> hasattr(T, 'volume')
False