Manifold: the main class

class snappy.Manifold

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, where each tetrahedron has 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 and some orbifolds. Here’s a quick example:

>>> M = Manifold('9_42')
>>> M.volume()
4.05686022424
>>> M.cusp_info('shape')
[(-4.27893631592+1.95728679750j)]

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}, [a_0, B_1])’ or ‘Splitting(S_{1,0}, [b_1, A_0], [a_0,B_1])’

    Type “twister?” for more details.

  • 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 string containing the contents of a SnapPea triangulation or link projection file.

canonize

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

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

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

chern_simons

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

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

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

>>> M.chern_simons().accuracy
9

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.1532041333
>>> M.dehn_fill( (1,2) )
>>> M.chern_simons()
0.077317871386

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

complex_volume
Returns the complex volume, i.e.
volume + i 2 pi^2 (chern simons)
>>> M = Manifold('5_2')
>>> M.complex_volume()
(2.8281220883-3.0241283765j)
copy

Returns a copy of the manifold

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

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 = snappy.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()
7.9999999999999947

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

sage: N3 = M.cover(f.Kernel())
sage: N3.volume()/M.volume()
167.99999999999858

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())
sage: Q, f = G.pQuotient(5, 1, nvals = 2)
sage: M.cover(f.Kernel()).volume()
10.149416064096533
sage: h = G.SimpleQuotients(1, 11, 2, 10^4)[1,1]
sage: N4 = M.cover(h)
sage: N2 == N4
True
covers

M.covers(degree, method=None)

Returns a list of Manifolds corresponding to all of the finite covers of the given degree.

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

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

You can also look just at cyclic covers, which is much faster.

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

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

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 coeffients (M, L) = (1.0, 2.0)

To get more detailed information about the cusp, we do

>>> c = M.cusp_info(0)
>>> c.shape
(0.110445017621+0.946770978498j)
>>> c.modulus
(-0.12155871955249957+1.042041282932261j)
>>> c.keys()
['index', 'holonomies', 'holonomy_accuracy', 'shape', 'filling', 'shape_accuracy', 'modulus', 'is_complete', '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.598830888594+1.098125481710j), (0.89824633289+1.49440443102j))

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.110445017621+0.946770978498j),
 Cusp 1 : torus cusp with Dehn filling coeffients (M, L) = (1.0, 2.0),
 Cusp 2 : torus cusp with Dehn filling coeffients (M, L) = (3.0, 2.0)]
>>> M.cusp_info('is_complete')
[True, False, False]
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.3247595264191645
>>> len(CN.horoballs(0.01))
178
>>> CN.view()  # Opens 3-d picture of the horoballs 
dehn_fill

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
    L408001(0,0)(0,0)(2,3)(0,0)
    
  • Fill the last cusp:

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

    >>> M.dehn_fill( [ (3,0), (1, -4) ])
    >>> M
    L408001(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

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-graphic of the DirichletDomain.  

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
drill

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 chosing which one to drill out.

>>> M = Manifold('v3000')
>>> N = M.drill(0, max_segments=3)
>>> (M.num_cusps(), N.num_cusps())
(1, 2)
dual_curves

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.562399148646-2.81543088521j),
   1: orientation-preserving curve of length (1.12479829729+0.652323536768j),
   2: orientation-preserving curve of length (1.26080401747+1.97804689023j),
   3: orientation-preserving curve of length (1.58826932598+1.67347167369j),
   4: orientation-preserving curve of length (1.68719744594+2.81543088521j)]

Each curve is returned as an info object with these keys

>>> curves[0].keys()
['index', 'filled_length', 'complete_length', '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.562399148646-2.81543088521j)]
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}
filled_triangulation

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
>>> type(N) == Triangulation
True

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

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')
matrix([[-2.50000000+2.59807621j,  6.06217783+0.5j       ],
        [-0.86602540+2.5j       ,  4.00000000-1.73205081j]])

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

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) = 1 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()
matrix([[ 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)]
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_isometric_to

Returns True if M and N are isometric, False otherwise.

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

length_spectrum

M.length_spectrum(cutoff=1.0)

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

name

Return the name of the triangulation.

>>> M = Triangulation('4_1')
>>> M.name()
'L104001'
num_cusps

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 manifold, returns the 2-fold cover which is orientable.

>>> X = Manifold('x123')
>>> Y = X.orientation_cover()
>>> (X.is_orientable(), Y.is_orientable())
(False, True)
>>> Y
x123~(0,0)(0,0)

Brings the corresponding link editor window to the front, if there is one.

randomize

Perform random Pachner moves on the underlying triangulation.

>>> 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()
>>> cs + M.chern_simons()
0.0
save

Save the triangulation as a SnapPea triangulation file.

>>> M = Triangulation('m004')
>>> M.save('fig-eight.tri')
set_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

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

This method changes these markings.

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

    >>> M = Manifold('5_2')
    >>> M.cusp_info('shape')
    [(-2.49024466751+2.97944706648j)]
    >>> M.set_peripheral_curves('shortest')
    >>> M.cusp_info('shape')
    [(-0.49024466751+2.97944706648j)]
    

    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

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

M.set_tetrahedra_shapes(shapes, fillings=[(1,0)]):

Replaces the tetrahedron shapes with those in the list ‘shapes’ and sets the Dehn filling coefficients as specified by the fillings argument.

simplify

Try to simplify the triangulation by doing Pachner moves.

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

Returns the type of the current solution to the gluing equations, basically a summary of how degenerate the solution is. The possible answers are:

  • ‘not attempted’
  • ‘all tetrahedra positively oriented’ aka ‘geometric_solution’ Should correspond to a genuine hyperbolic structure
  • ‘contains negatively oriented tetrahedra’ aka ‘nongeometric_solution’ Probably correponds to a hyperbolic structure but some simplices have reversed orientiations.
  • ‘contains flat tetrahedra’ Contains some tetrahedra with shapes in R - {0, 1}.
  • ‘contains degenerate tetrahedra’ Some shapes are close to {0,1, or infinity}.
  • ‘unrecognized solution type’
  • ‘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'
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

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

tetrahedra_shapes

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 choosen so as to normalize the triangle.

>>> M = Manifold('m015')
>>> M.tetrahedra_shapes(part='rect')
[(0.6623589786223731+0.5622795120623011j), (0.662358978622373+0.5622795120623011j), (0.6623589786223729+0.562279512062301j)]
>>> M.tetrahedra_shapes()
[{'accuracies': (11, 11, 12, 11), 'log': (-0.14059978716148094+0.7038577213014763j), 'rect': (0.6623589786223731+0.5622795120623011j)},
 {'accuracies': (11, 11, 11, 11), 'log': (-0.14059978716148103+0.7038577213014764j), 'rect': (0.662358978622373+0.5622795120623011j)},
 {'accuracies': (11, 11, 11, 11), 'log': (-0.14059978716148125+0.7038577213014764j), 'rect': (0.6623589786223729+0.562279512062301j)}]
volume

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

>>> M.volume().accuracy
10
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.02988321282
without_hyperbolic_structure

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

>>> M = Manifold('9_42')
>>> T = M.without_hyperbolic_structure()
>>> type(T)
<type 'snappy.SnapPy.Triangulation'>

Previous topic

The snappy module and its classes

Next topic

Triangulation

This Page