Additional Classes

SnapPy is a Cython wrapping of Jeff Weeks’ SnapPea kernel.

The module defines the following classes:

Triangulation, Manifold, ManifoldHP, AbelianGroup, FundamentalGroup, HolonomyGroup, HolonomyGroupHP, DirichletDomain, DirichletDomainHP, CuspNeighborhood, CuspNeighborhoodHP, SymmetryGroup, AlternatingKnotExteriors, NonalternatingKnotExteriors, SnapPeaFatalError, InsufficientPrecisionError, pari, twister, OrientableCuspedCensus, NonorientableCuspedCensus, OrientableClosedCensus, NonorientableClosedCensus, LinkExteriors, CensusKnots, HTLinkExteriors, TetrahedralOrientableCuspedCensus, TetrahedralNonorientableCuspedCensus, OctahedralOrientableCuspedCensus, OctahedralNonorientableCuspedCensus, CubicalOrientableCuspedCensus, CubicalNonorientableCuspedCensus, DodecahedralOrientableCuspedCensus, DodecahedralNonorientableCuspedCensus, IcosahedralNonorientableClosedCensus, IcosahedralOrientableClosedCensus, CubicalNonorientableClosedCensus, CubicalOrientableClosedCensus, DodecahedralNonorientableClosedCensus, DodecahedralOrientableClosedCensus, Crossing, Strand, Link, Tangle, RationalTangle, ZeroTangle, InfinityTangle, IdentityBraid, random_link, DTcodec.

AbelianGroup

class snappy.AbelianGroup(presentation=None, elementary_divisors=[])

An AbelianGroup object represents a finitely generated abelian group, usually the first homology group of a snappy Manifold.

Instantiate an abelian group by its elementary divisors:

>>> A = AbelianGroup(elementary_divisors=[5,15,0,0])
>>> A
Z/5 + Z/15 + Z + Z
>>> A[0]
5

Alternatively, instantiate an abelian group as AbelianGroup(P) where P is a presentation matrix given as a list of lists of integers. Snappy stores an abelian group as a list of elementary divisors:

>>> B = AbelianGroup([[1,3,2],[2,0,6]])
>>> B
Z/2 + Z
>>> B.elementary_divisors()
[2, 0]
>>> B[1]
0
>>> len(B)
2
>>> B.rank()
2
>>> B.betti_number()
1
>>> B.order()
'infinite'
betti_number()

The rank of the maximal free abelian subgroup.

elementary_divisors()

The elementary divisors of this finitely generated abelian group.

order()

The order of the group. Returns the string ‘infinite’ if the group is infinite.

rank()

The rank of the group.

FundamentalGroup

class snappy.HolonomyGroup

A FundamentalGroup represents a presentation of the fundamental group of a SnapPea Triangulation. Group elements are described as words in the generators a,b,…, where the inverse of a is denoted A. Words are represented by python strings (and the concatenation operator is named ‘+’, according to Python conventions).

Instantiate via T.fundamental_group(), where T is a triangulation:

>>> T = Triangulation('m125')
>>> T.fundamental_group()
Generators:
   a,b
Relators:
   aabaBBAABAbb
>>> type(T.fundamental_group()) 
<class 'SnapPy.FundamentalGroup'>

A HolonomyGroup is a FundamentalGroup with added structure consisting of a holonomy representation into O(3,1), and an arbitrarily chosen lift of the holonomy representation to SL(2,C). The holonomy is determined by the shapes of the tetrahedra, so a HolonomyGroup is associated to a Manifold, while a Triangulation only has a FundamentalGroup:

Instantiate via M.fundamental_group(), where M is a Manifold:

>>> M = Manifold('m125')
>>> G = M.fundamental_group()
>>> G
Generators:
   a,b
Relators:
   aabaBBAABAbb
>>> type(G) 
<class 'SnapPy.HolonomyGroup'>

In the class HolonomyGroup, methods are provided to evaluate the representations on a group element. Other methods are shared with the FundamentalGroup class.

>>> G.O31('a') 
[    2.50000000000000   -0.500000000000000    -2.12132034355964   -0.707106781186547]
[   0.500000000000002   -0.500000000000001   -0.707106781186549    0.707106781186547]
[  -0.707106781186548   -0.707106781186547     1.00000000000000                    0]
[    2.12132034355964   -0.707106781186548    -2.00000000000000    -1.00000000000000]
>>> G.SL2C('aaB') 
[-1.00000000000000 + 4.00000000000000*I 2.12132034355964 - 0.707106781186545*I]
[2.12132034355964 - 0.707106781186549*I -1.00000000000000 - 1.00000000000000*I]
>>> G.complex_length('ab') 
1.06127506190504 - 2.23703575928741*I
O31(word)

Return the image of the element represented by the input word under the holonomy representation, where Isom(H^3) is identified with SO(3,1).

SL2C(word)

Return the image of the element represented by the input word under some SL(2,C) representation that lifts the holonomy representation. Note: the choice of lift is not guaranteed to vary continuously when filling coefficients are changed.

character_variety_vars_and_polys(as_ideal=False)

Returns a list of variables and a list polynomials where the polynomials generate the ideal defining the SL(2, C) character variety of this group. Each variable is of the form “Tw” where “w” is a word in the generators and “Tw” represents the trace function of that word.

>>> H = Manifold('dLQacccbjkg')  # Hopf link exterior.
>>> G = H.fundamental_group()
>>> vars, polys = G.character_variety_vars_and_polys()
>>> vars
[Ta, Tb, Tab]
>>> polys    
[Ta^3 - Tab*Tb*Ta^2 + (Tb^2 + (Tab^2 - 4))*Ta,
 Ta^2 - Tab*Tb*Ta + (Tb^2 + (Tab^2 - 4))]

When used inside Sage, you can ask for the answer as a proper ideal:

sage: M = Manifold('m000')  # Gieseking manifold
sage: G = M.fundamental_group()
sage: I = G.character_variety_vars_and_polys(as_ideal=True)
sage: I
Ideal (-Ta^3*Tb^2*Tab + Ta^4*Tb + Ta^2*Tb^3 + Ta^2*Tb*Tab^2 + Ta*Tb^2*Tab - 5*Ta^2*Tb - Tb^3 - Tb*Tab^2 + Ta*Tab - Ta + 3*Tb, Tb*Tab - Ta - Tb, -Ta^2*Tb^2*Tab + Ta^3*Tb + Ta*Tb^3 + Ta*Tb*Tab^2 - 4*Ta*Tb + Tab - 2) of Multivariate Polynomial Ring in Ta, Tb, Tab over Rational Field
sage: I.dimension()
1
complex_length(word)

Return the complex length of the isometry represented by the input word.

gap_string()

Returns a string which will define this group within GAP:

>>> M = Manifold('b++LLR')
>>> G = M.fundamental_group()
>>> G
Generators:
   a,b
Relators:
   aaaaBAbbAB
>>> G.gap_string()
'CallFuncList(function() local F, a, b; F := FreeGroup("a", "b"); a := F.1; b := F.2;   return F/[a*a*a*a*b^-1*a^-1*b*b*a^-1*b^-1]; end,[])'
>>> G.magma_string()
'Group<a,b|a*a*a*a*b^-1*a^-1*b*b*a^-1*b^-1>'
generators()

Return the letters representing the generators in the presentation.

>>> M = Manifold('9_42')
>>> G = M.fundamental_group() #Presentation simplified by default
>>> G
Generators:
   a,b
Relators:
   aaaabbABBBAbb
>>> H = M.fundamental_group(False) #Unsimplified presentation
>>> H
Generators:
   a,b,c,d,e
Relators:
   ECbC
   dEb
   dAcaB
   dbaE

SnapPy stores a FundamentalGroup as a presentation of the group. The following commands demonstrate how generators in the unsimplified and simplified presentations above correspond:

>>> G.generators()
['a', 'b']
>>> H.generators()
['a', 'b', 'c', 'd', 'e']
>>> G.generators_in_originals()
['BABcBcbCABcBcbCCbCba', 'BcBCbCbab']
>>> G.original_generators()
['BBAbba', 'A', 'AB', 'abba', 'bba']
>>> G.num_generators()
2
>>> G.num_original_generators()
5
>>> H.num_generators()
5
generators_in_originals(verbose_form=False, raw_form=False, as_int_list=False)

Return the current generators in terms of the original geometric generators. Note that by default fundamental_group() returns a simplified presentation of the group.

If the flag “raw_form” is set to True, it returns a sequence of instructions for expressing the current generators in terms of the original ones. This is sometimes much more concise, though the format is somewhat obscure. See the source code of this function in fundamental_group.pyx for details.

>>> M = Manifold('K7_1')
>>> G = M.fundamental_group()
>>> G.generators_in_originals()   
['DBcACbABcaCbDBcACbABcaCbACbaBcaCbd',
 'DBcACbABcaCbDBcACbABcaCbCbDBcACbABcaBcACbaBcaCbdACbaBcaCbdBcBcACbaBcaCbd']
longitude(int which_cusp=0, as_int_list=False)

Returns a word representing a conjugate of the current longitude for the given cusp. Guaranteed to commute with the meridian for the same cusp. Note: for Klein bottle cusps, the longitude must be defined carefully.

>>> G = Manifold('m004').fundamental_group()
>>> G.longitude(0)
'aBAbABab'
>>> G.longitude()   # shortcut for the above.
'aBAbABab'
magma_string()

Returns a string which will define this group within MAGMA.

meridian(int which_cusp=0, as_int_list=False)

Returns a word representing a conjugate of the current meridian for the given cusp. Guaranteed to commute with the longitude for the same cusp.

>>> G = Manifold('m125').fundamental_group()
>>> G.meridian(0)
'aaba'
>>> G.meridian(-1)  # The last cusp
'baaba'
num_generators()

Return the number of generators for the presentation.

num_original_generators()

Return the number of geometric generators (before simplification).

num_relators()

Return the number of generators for the presentation.

original_generators(verbose_form=False, as_int_list=False)

Return the original geometric generators (before simplification) in terms of the current generators.

>>> M = Manifold('v0000')
>>> G = M.fundamental_group()
>>> G.original_generators(as_int_list=True)
[[1], [-1, -2, 1, 2], [-1, 2], [-2, -1, 2], [2]]
peripheral_curves(as_int_list=False)

Returns a list of meridian-longitude pairs for all cusps.

>>> G = Manifold('m125').fundamental_group()
>>> G.peripheral_curves()
[('aaba', 'abb'), ('baaba', 'Ba')]
relators(verbose_form=False, as_int_list=False)

Return a list of words representing the relators in the presentation.

If the optional argument verbose_form is True, then the relator is returned in the form “a*b*a^-1*b^-1” instead of “abAB”.

sage()

Returns the corresponding Sage FinitelyPresentedGroup

SymmetryGroup

class snappy.SymmetryGroup

A SymmetryGroup is a group of self-isometries of hyperbolic 3-manifold. Instantiate as follows:

>>> M = Manifold('m004')
>>> M.symmetry_group()
D4
abelian_description()

If the symmetry group is abelian, return it as an AbelianGroup

>>> S = Manifold('v3379').symmetry_group()
>>> S.abelian_description()
Z/2 + Z/2 + Z/2
abelianization()

Return the abelianization of the symmetry group

>>> S = Manifold('m004').symmetry_group()
>>> S.abelianization()
Z/2 + Z/2
center()

Return the center of the symmetry group

>>> S = Manifold('m004').symmetry_group()
>>> S.center()
Z/2
commutator_subgroup()

Return the commutator subgroup of the SymmetryGroup

>>> S = Manifold('m004').symmetry_group()
>>> S
D4
>>> S.commutator_subgroup()
Z/2
direct_product_description()

If the SymmetryGroup is a nontrivial direct product with at least one nonabelian factor, return a pair of SymmetryGroups consisting of the (two) factors.

>>> S = Manifold('s960').symmetry_group()
>>> S.direct_product_description()
(Z/4, D3)
is_S5()

Returns whether the group is the symmetric group on five things.

is_abelian()

Return whether the symmetry group is abelian.

>>> S = Manifold('m004').symmetry_group()
>>> S.is_abelian()
False
is_amphicheiral()

Return whether the manifold has an orientation reversing symmetry.

>>> S = Manifold('m004').symmetry_group()
>>> S.is_amphicheiral()
True
is_dihedral()

Return whether the symmetry group is dihedral.

>>> S = Manifold('m004').symmetry_group()
>>> S.is_dihedral()
True
is_direct_product()

Return whether the SymmetryGroup is a nontrivial direct product with at least one nonabelian factor.

>>> S = Manifold('s960').symmetry_group()
>>> S.is_direct_product()
True
>>> S
Z/4 x D3
is_full_group()

Return whether the full symmetry group has been found.

>>> S = Manifold('m004').symmetry_group()
>>> S.is_full_group()
True
is_invertible_knot()

Return whether a one-cusped has a symmetry that acts on the cusp via the matrix -I.

>>> S = Manifold('m015').symmetry_group()
>>> S.is_invertible_knot()
True
is_polyhedral()

Returns whether the symmetry group is a (possibly binary) polyhedral group.

isometries()

Return a detailed list of all the isometries in the symmetry group.

>>> S = Manifold('s959').symmetry_group()
>>> isoms = S.isometries()
>>> isoms[8]
0 -> 1   1 -> 0
[-1 -1]  [ 0  1]
[ 1  0]  [-1 -1]
Does not extend to link
multiply_elements(i, j)

Returns the product of group elements i and j. The convention is that products of symmetries read right to left. That is, the composition (symmetry[i] o symmetry[j]) acts by first doing symmetry[j], then symmetry[i].

>>> S = Manifold('m004').symmetry_group()
>>> S.multiply_elements(2, 3)
1
order()

Return the order of the symmetry group

>>> S = Manifold('s000').symmetry_group()
>>> S.order()
4
polyhedral_description()

If the symmetry group is a (possibly binary) polyhedral group, return a description of it.

DirichletDomain

class snappy.DirichletDomain

A DirichletDomain object represents a Dirichlet Domain of a 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.

Instantiate as M.dirichlet_domain() where M is a Manifold to obtain a Dirichlet Domain centered at a point which maximizes injectivity radius.

Other options can be provided to customize the computation, with the default values shown here

>>> M = Manifold('m003(3,-4)')
>>> M.dirichlet_domain(vertex_epsilon=10.0**-8, displacement = [0.0, 0.0, 0.0], centroid_at_origin=True, maximize_injectivity_radius=True)
40 finite vertices, 0 ideal vertices; 60 edges; 22 faces

You can also create a Dirichlet Domain from a file listing matrix generators for the group, in SnapPea’s “% Generator” format, via

D = DirichletDomain(generator_file=’test.gens’)

Or from matrices:

>>> G = Manifold('m004').fundamental_group(False)
>>> matrices = [ G.O31('a'), G.O31('b'), G.O31('c') ] # Note: some of the matrices contain (near) 0 entries and thus this tests that Object2Real converts small numbers fromatted by pari as "1.0 E-10" (note the pace before "E") correctly when not in SageMath.
>>> DirichletDomain(O31_generators = matrices,
...                 maximize_injectivity_radius = False)
8 finite vertices, 2 ideal vertices; 20 edges; 12 faces

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

edge_list()

Return a list of edges, each represented as a dictionar with keys ‘tail_vertex_index’, ‘tip_vertex_index’, ‘edge_class’.

The index (into vertex_data_list()) to the two vertices at the end of the edge are stored in ‘tail_vertex_index’ and ‘tip_vertex_index’. The index of the edge class this edge belongs to is stored in ‘edge_class’.

export_stl(filename, model='klein', cutout=False, num_subdivisions=3, shrink_factor=0.9, cutoff_radius=0.9, callback=None)

Export the Dirichlet domain as an stl file suitable for 3d printing.

Arguments can be given to modify the model produced:

  • model=’klein’ - (alt. ‘poincare’) the model of HH^3 to use.

  • cutout=False - remove the interior of each face

  • shrink_factor=0.9 - the fraction to cut out of each face

  • cuttoff_radius=0.9 - maximum rescaling for projection into Poincare model

  • num_subdivision=3 - number of times to subdivide for the Poincare model

For printing domains in the Poincare model, cutoff_radius is critical for avoiding infinitely thin cusps, which cannot be printed.

This can take a long time for finely subdivided domains. So we call UI_callback every so often if it is not None.

>>> D = Manifold('m004').dirichlet_domain()
>>> D.export_stl('fig-eight-klein.stl')     
>>> D.export_stl('fig-eight-poincare.stl', model='poincare')     
>>> D.export_stl('fig-eight-klein-wireframe.stl', cutout=True)     
>>> D.export_stl('fig-eight-poincare-wireframe.stl', model='poincare', cutout=True)     
face_list()

Return a list of faces, each represented as a dictionary with keys ‘vertices’, ‘distance’, ‘closest’, ‘hue’, ‘vertex_indices’, ‘edge_indices’, ‘vertex_image_indices’, ‘edge_image_indices’, ‘edge_orientations’.

The distance from the origin is the value for ‘distance’, and the value for ‘closest’ is the orthogonal projection of the origin to the plane containing the face. The vertices of each face are listed in clockwise order, as viewed from outside the polyhedron.

The coordinates of vertices are stored in ‘vertices’ and the corresponding index into vertex_data_list() is stored in ‘vertex_index’. The indices (in edge_list()) to the edges of the face (also in clockwise order) are stored in ‘edge_indices’ such that the first edge is adjacent to the first and second vertex. The respective value in ‘edge_orientations’ is +/-1 to indicate whether the orientation of the edge induced from the orientation of the face is the same or opposite than the edges orientation.

To find the image of a vertex or edge adjacent to a face under the pairing matrix for this face, lookup the index in ‘vertex_image_indices’, respectively, ‘edge_image_indices’ at the respective position.

in_radius()

Return the radius of the largest inscribed sphere.

length_spectrum_dicts(cutoff_length=1.0, full_rigor=True, multiplicities=True, user_radius=0.0, grouped=True)

Return a list of info objects describing the short geodesics up to the specified cutoff length. The keys are ‘length’, ‘parity’, ‘topology’, and ‘multiplicity’. The length is the complex length; the parity specifies whether orientation is preserved; and topology distinguishes between circles and mirrored intervals. Finally, the key ‘matrix’ in the fundamental group realizing this element.

>>> M = Manifold('m004(1,2)')
>>> D = M.dirichlet_domain(maximize_injectivity_radius=False)
>>> lengths = D.length_spectrum_dicts()
>>> len(lengths)
2
>>> lengths[0].matrix in D.pairing_matrices()
True

If the flag ‘grouped’ is False, then each geodesic is returned as a separate item rather than collating by (length, parity, topology). If the flag ‘multiplicities’ is False, then the geodesics are collated but the multiplicity of each item is set to 0.

>>> M = Manifold('m003(-3, 1)')
>>> D = M.dirichlet_domain()
>>> [g.multiplicity for g in D.length_spectrum_dicts()]
[3, 3]
>>> [g.multiplicity for g in D.length_spectrum_dicts(grouped=False)]
[1, 1, 1, 1, 1, 1]
manifold()

Returns a Manifold computed directly from the Dirichlet domain, regarded as polyhedron with faces identified in pairs. Only works if this gives a manifold not an orbifold.

>>> M = Manifold('7_3')
>>> D = M.dirichlet_domain()
>>> A = D.manifold()
>>> M.is_isometric_to(A)
True
num_edges()

Return the number of edges.

num_faces()

Return the number of faces.

num_finite_vertices()

Return the number of finite (non-ideal) vertices.

num_ideal_vertices()

Return the number of ideal vertices.

num_vertices()

Return the number of vertices.

out_radius()

Return the radius of the smallest circubscribed sphere.

pairing_matrices()

Returns a list of the O31Matrices which pair the faces of this DirichletDomain.

>>> M = Manifold('s345')
>>> D = M.dirichlet_domain()
>>> matrices = D.pairing_matrices()
>>> D1 = DirichletDomain(O31_generators=matrices)
>>> N = D1.manifold()
>>> M.is_isometric_to(N)
True
pairing_words()

Group elements which pair the faces of this DirichletDomain as words in the original generators:

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

Requires that DirichletDomain was computed with include_words = True.

save(filename)

Save the Dirichlet domain as a text file in “% Generators” format.

>>> from snappy.number import Number
>>> acc, Number._accuracy_for_testing = Number._accuracy_for_testing, None
>>> M = Manifold('m125')
>>> D = M.dirichlet_domain()
>>> from tempfile import NamedTemporaryFile
>>> f = NamedTemporaryFile(delete=False)
>>> f.close()
>>> D.save(f.name)
>>> E = DirichletDomain(generator_file=f.name); E
30 finite vertices, 2 ideal vertices; 50 edges; 20 faces
>>> os.unlink(f.name)
>>> from pickle import dumps, loads
>>> E = loads(dumps(D)); E
30 finite vertices, 2 ideal vertices; 50 edges; 20 faces
>>> Number._accuracy_for_testing = acc
spine_radius()

Return the infimum of the radii (measured from the origin) of all spines dual to the Dirichlet domain.

triangulation()

Returns a Triangulation computed directly from the Dirichlet domain, regarded as polyhedron with faces identified in pairs. Only works if this gives a manifold not an orbifold.

>>> M = Manifold('7_3')
>>> D = M.dirichlet_domain()
>>> B = D.triangulation()
>>> M.is_isometric_to(B.with_hyperbolic_structure())
True
vertex_list(details=False)

Return a list of the coordinates of the vertices. These are the three space coordinates of a point in the time=1 slice of Minkowski space. That is to say, these are the coordinates of the image of the point under projection into the Klein model.

If details = True is passed, returns a list of vertices, each represented by a dictionary with keys ‘position’, ‘ideal’, ‘vertex_class’. The coordinates are the value for ‘position’. The index of the vertex class this vertex belongs to is the value for ‘vertex_class’. The value for ‘ideal’ is True if the vertex is an ideal point.

view()
volume()

Returns the approximate volume of the DirichletDomain. Because matrices in O(3,1) tend to accumulate roundoff error, it’s hard to get a good bound on the accuracy of the computed volume. Nevertheless, the kernel computes the best value it can, with the hope that it will aid the user in recognizing manifolds defined by a set of generators.

CuspNeighborhood

class snappy.CuspNeighborhood

A CuspNeighborhood object represents an equivariant collection of disjoint horoballs that project to cusp neighborhoods.

Instantiate as M.cusp_neighborhood()

Ford_domain(which_cusp=0, high_precision=False)

Return a list of pairs of complex numbers describing the endpoints of the segments obtained by projecting the edges of the Ford domain to the xy-plane in the upper half space model.

If the high_precision flag is set to False (the default), the coordinates are Python complex numbers. Otherwise they are SnapPy Numbers.

all_translations(verified=False, bits_prec=None)

Returns the (complex) Euclidean translations of the meridian and longitude for each cusp measured with respect to the cusp neighborhood.

The result is a list of pairs, the second entry corresponding to a longitude is always real:

>>> M = Manifold("v3227")
>>> N = M.cusp_neighborhood()
>>> N.all_translations() 
[(-0.152977162509284 + 0.747697694854404*I, 0.868692062725708), (-0.152977162509284 + 0.747697694854404*I, 0.868692062725708), (0.0961611977895952 + 0.725536253181650*I, 0.895226186134782)]

Often, one is interested in making the cusp neighborhoods as large as possible first:

>>> N.set_displacement(100,0)
>>> N.set_displacement(100,1)
>>> N.set_displacement(100,2)
>>> N.all_translations() 
[(-0.477656250512815 + 2.33461303362557*I, 2.71240613125259), (-0.259696455247511 + 1.26930345526993*I, 1.47470541152065), (0.131389112265699 + 0.991330873713731*I, 1.22318540718077)]

This can also be achieved by Manifold.cusp_translations() which would have made a different choice of disjoint cusp neighborhoods though:

>>> M.cusp_translations() 
[(-0.315973594129651 + 1.54436599614183*I, 1.79427928161946), (-0.315973594129649 + 1.54436599614182*I, 1.79427928161946), (0.198620491993677 + 1.49859164484929*I, 1.84908538602825)]

This method supports arbitrary precision

>>> from snappy.number import Number
>>> N.set_displacement(1.125, 0)
>>> N.set_displacement(0.515625, 1)
>>> N.set_displacement(0.3125, 2)
>>> N.all_translations(bits_prec = 120) 
[(-0.47120283346076781167174343474008914 + 2.3030710375877078211095122873223488*I, 2.6757599281290843845710310925394911), (-0.25618853688042434043044508297577899 + 1.2521580040549576537090841783446072*I, 1.4547854392045669515377748986943560), (0.13143677360753666862808198126761923 + 0.99169047854575721271560179767750893*I, 1.2236291171413362101960100623801910)]

and can return verified intervals

sage: N.all_translations(verified = True) # doctest: +NUMERIC9
[(-0.47120283346? + 2.30307103759?*I, 2.67575992813?), (-0.256188536881? + 1.252158004055?*I, 1.454785439205?), (0.131436773608? + 0.991690478546?*I, 1.2236291171413?)]
sage: N.all_translations(verified = True, bits_prec = 120) # doctest: +NUMERIC30
[(-0.4712028334607678116717434347401? + 2.3030710375877078211095122873224?*I, 2.6757599281290843845710310925395?), (-0.25618853688042434043044508297578? + 1.25215800405495765370908417834461?*I, 1.454785439204566951537774898694356?), (0.131436773607536668628081981267619? + 0.991690478545757212715601797677509?*I, 1.223629117141336210196010062380191?)]

that are guaranteed to contain the true translations of disjoint cusp neighborhoods (the element corresponding to a longitude is always in a RealIntervalField). The verified translations might correspond to cusp neighborhoods smaller than the given ones to be able to verify that they are disjoint.

Remark: Since the code is (potentially) non-deterministic, the result of

[ N.all_translations(verified = True)[i] for i in range(M.num_cusps()) ]

is not verified to correspond to disjoint cusp neighborhoods.

check_index(which_cusp)

Raises an IndexError if the cusp index is invalid.

get_displacement(which_cusp=0)

Return the displacement of the horospherical boundary of the specified cusp. The displacement is the hyperbolic distance that the horospherical boundary has been displaced from its “home” position, at which the area of the boundary is 3sqrt(3)/8. (The translates of all of the horospheres are guaranteed to be pairwise disjoint when each cusp has displacement 0.)

get_tie(which_cusp)

Return True if the specified cusp is a member of the tied group. The displacements of the tied cusps are all the same.

horoballs(cutoff=0.1, which_cusp=0, full_list=True, high_precision=False)

Return a list of dictionaries describing the horoballs with height at least cutoff. The keys are ‘center’, ‘radius’, ‘index’.

If the high_precision flag is set to the default value False, these are Python complexes and floats. Otherwise they are SnapPy Numbers.

manifold()

Return a Manifold built from the current canonical triangulation.

max_reach()

Return the maximum reach over all cusps.

num_cusps()

Return the number of cusps.

original_index(which_cusp)

Returns the index by which the Manifold identifies this cusp.

reach(which_cusp=0)

Return the displacement at which the specified cusp neighborhood bumps into itself. (This is twice the distance between nearest horoball lifts.)

set_displacement(new_displacement, which_cusp=0)

Set the displacement of the specified cusp.

set_tie(which_cusp, new_tie)

Mark the specified cusp as a member of the tied group.

stopper(which_cusp)

Return the index of the cusp which will be the first one that the specified cusp neighborhood bumps into. (Assumes the other displacements are fixed.)

stopping_displacement(which_cusp=0)

Return the displacement at which the specified cusp neighborhood bumps into itself or another cusp neighborhood. (Assumes the other displacements are fixed.)

topology(which_cusp=0)

Return the topological type of the specified cusp.

translations(which_cusp=0)

Return the (complex) Euclidean translations of the meridian and longitude of the specified cusp.

Also see CuspNeighborhood.all_translations() which supports arbitrary precision and verified results.

triangulation(which_cusp=0, high_precision=False)

Return a list of dictionaries describing the endpoints of the segments obtained by projecting the edges of the triangulation dual to the Ford domain into the xy-plane in the upper half space model. The keys are ‘endpoints’ and ‘indices’.

view(which_cusp=0, cutoff=None)

Create a 3D picture of the horoball packing. One can specify which cusp to put at infinity and how large of horoballs to look at, e.g.

>>> M = Manifold('m125')
>>> C = M.cusp_neighborhood()
>>> C.view(which_cusp = 1, cutoff=0.2)   
volume(which_cusp=0)

Return the volume of the horoball neighborhood of the specified cusp.