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 and Klein-bottle 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
. SeeManifold.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: bool = False, interval_bits_precs: Sequence[int] = [53, 212], exact_bits_prec_and_degrees: Sequence[Tuple[int, int]] = [(212, 10), (1000, 20), (2000, 20)], verbose: bool = False) Triangulation | Manifold
Returns a triangulation canonically associated to the hyperbolic manifold. That is, the triangulation is (up to combinatorial isomorphism relabeling the tetrahedra and vertices) completely determined by the isometry type of the hyperbolic manifold.
Manifolds with incomplete cusps are rejected (unlike in the case of
isometry_signature
).We now describe the canonical retriangulation. If all cells of the canonical cell decomposition (defined by Epstein and Penner ‘88) are tetrahedral,
canonical_retriangulation
simply returns that ideal triangulation as aManifold
. Here is an example:>>> M = Manifold("m015") >>> K = M.canonical_retriangulation() >>> K.has_finite_vertices() False >>> K.solution_type() 'all tetrahedra positively oriented'
If there are non-tetrahedral cells,
canonical_retriangulation
subdivides the canonical cell decomposition. It introduces a finite vertex for each canonical cell resulting in aTriangulation
. Here is an example where the canonical cell is a cube:>>> M = Manifold("m412") >>> K = M.canonical_retriangulation() >>> K.has_finite_vertices() True
The canonical retriangulation can be used to find the symmetries of a single manifold. It also can compute the isometries between two manifolds. We do this using
isomorphisms_to
:>>> M = Manifold("5_2").canonical_retriangulation() >>> N = Manifold("m015").canonical_retriangulation() >>> M.isomorphisms_to(M) [0 -> 0 [1 0] [0 1] ... >>> M.isomorphisms_to(N) [0 -> 0 [-1 2] [ 0 -1] ...
The canonical retriangulation is also the basis for the
isometry_signature
.Subdivision
If the canonical cell decomposition has a non-tetrahedral cell, the method subdivides. You can think of the subdivision in either of the following (equivalent) ways:
A coarsening of the barycentric subdivision with only a quarter of the number of tetrahedra. That is, take the barycentric subdivision and merge the four tetrahedra adjacent to a barycentric edge connecting an edge midpoint to a face midpoint.
Taking the double suspension of each face (which is an ideal n-gon) about the centers of the two neighboring 3-cells. Then split each such topological “lens” into n tetrahedra along its central axis.
Verified computations
While the canonical retriangulation is combinatorial, some intermediate computations are numerical. Thus, if
verified = False
, floating-point issues can arise. (Arguably this gave rise to a mistake in the non-orientable census.x101
andx103
were later identified as the same by Burton ‘14.)The method can be made verified by passing
verified = True
:sage: M = Manifold("v2986") sage: K = M.canonical_retriangulation(verified = True) sage: K.has_finite_vertices() # Cell decomposition verified to be tetrahedral False sage: K.triangulation_isosig(decorated=False) # Verified isometry signature. 'jvLALQQdeefgihihiokcmmwwswg' sage: len(K.isomorphisms_to(K)) # Verified to have no (non-trivial) symmetries. 1
Interval arithmetic can only be used to verify the canonical cell decomposition if all cells are tetrahedral. For non-tetrahedral cells, the method automatically switches to exact methods to verify the canonical cell decomposition. That is, it uses snap-like methods (LLL-algorithm) to guess a representation of the shapes in the shape field. It then uses exact arithmetic to verify the shapes form a valid geometric structure and compute the necessary tilts to verify the canonical cell decomposition. Note that this can take a long time!
Here is an example where exact methods are used:
sage: M = Manifold("m412") sage: K = M.canonical_retriangulation(verified = True) sage: K.has_finite_vertices() # Has non-tetrahedral cell True
If the canonical retriangulation cannot be verified, an exception will be raised. (Note that this is new (and safer) in Version 3.2. Prior to that version,
Manifold.canonical_retriangulation()
could returnNone
instead.)Here is an example where we skip the (potentially lengthy) exact methods needed to verify a non-tetrahedral cell. The method fails (early and with an exception) since the cells are actually tetrahedral:
sage: M = Manifold("m412") sage: K = M.canonical_retriangulation(verified = True, exact_bits_prec_and_degrees = []) # doctest: +ELLIPSIS +IGNORE_EXCEPTION_DETAIL Traceback (most recent call last): ... snappy.verify.exceptions.TiltInequalityNumericalVerifyError: Numerical verification that tilt is negative has failed: ... < 0
- Parameters:
verified – Use verified computation.
interval_bits_precs – Only relevant if
verified = True
. A list of (increasing) precisions used to try to certify the canonical cell decomposition using intervals. Each precision is tried until we succeed. If none succeeded, we move on to exact methods.exact_bits_prec_and_degrees – Only relevant if
verified = True
. A list of pairs (precision, max degree) used when the LLL-algorithm is trying to find the defining polynomial of the shape field withListOfApproximateAlgebraicNumbers.find_field
. Each pair is tried until we succeed.verbose – Print information about the methods tried to compute and verify the canonical retriangulation.
- Returns:
If the canonical cell decomposition exists entirely of (hyperbolic ideal) tetrahedra, a
Manifold
with those tetrahedra. Otherwise, aTriangulation
that is a subdivision of the canonical cell decomposition.
- canonize() None
Change the triangulation to an arbitrary retriangulation of the canonical cell decomposition. See
canonical_retriangulation
to get the actual 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 actually not a retriangulation of the canonical cell decomposition.
- chern_simons(accuracy=False)
Returns the Chern-Simons invariant of the manifold (normalized by dividing it by \(2 \pi^2\)), 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
ValueError: The Chern-Simons invariant isn't 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 modulo \(i \pi^2\) which is given by
\[\text{vol} + i \text{CS}\]where \(\text{CS}\) is the (unnormalized) Chern-Simons invariant.
>>> M = Manifold('5_2') >>> M.complex_volume() 2.82812209 - 3.02412838*I
Note that
chern_simons
normalizes the Chern-Simons invariant by dividing it by \(2 \pi^2 = 19.7392...\)>>> M.chern_simons() -0.153204133297152
More examples:
>>> M.dehn_fill((1,2)) >>> M.complex_volume() 2.22671790 + 1.52619361*I >>> M = Manifold("3_1") # A non-hyperbolic example. >>> 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. Higher precision can be requested withbits_prec
:sage: M = Manifold("m015") sage: M.complex_volume(verified_modulo_2_torsion=True, bits_prec = 93) # doctest: +NUMERIC21 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) snappy.Manifold
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 listP
such such thatset(P) == set(range(d))
whered
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 SageMath, the permutations can also be of type
PermutationGroupElement
, in which case they act on the setrange(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: Optional[str] = None, cover_type: str = 'all') list[snappy.Manifold]
Returns a list of
Manifold
s 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: str = 'maximal', verified: bool = False, bits_prec: int | None = None)
Returns the maximal cusp area matrix \((A_{ij})\) where \(A_{ij}\) is defined as follows. Let \(C_i\) and \(C_j\) be the (open) cusp neighborhoods about cusp \(i\) and \(j\). Let \(A(C_i)\) and \(A(C_j)\) be the areas of \(C_i\) and \(C_j\), respectively. Then, \(C_i\) and \(C_j\) are embedded (if \(i = j\)) or disjoint (otherwise) if and only if \(A(C_i)A(C_j) \leq A_{ij}\).
Here is an example:
>>> M = Manifold("L6a5") >>> M.cusp_area_matrix() [27.9999999999996 7.00000000000000 7.00000000000000] [7.00000000000000 27.9999999999999 7.00000000000000] [7.00000000000000 7.00000000000000 28.0000000000001]
Faster lower bounds
This section can be skipped by most users!
Prior to SnapPy version 3.2, the algorithm to compute the maximal cusp area matrix was much slower and required
verified = True
and SageMath. Thus, in prior versions,method
defaulted totrigDependentTryCanonize
. This meant, that, by default,cusp_area_matrix()
only returned (some) lower bounds for the maximal cusp area matrix entries.These lower bounds can still be accessed:
>>> M.cusp_area_matrix(method = 'trigDependentTryCanonize') [21.4375000000000 7.00000000000000 7.00000000000000] [7.00000000000000 28.0000000000000 7.00000000000000] [7.00000000000000 7.00000000000000 28.0000000000000]
If
method = 'trigDependent'
ormethod = 'trigDependenyTryCanonize'
, the result is triangulation dependent or not even deterministic, respectively. Furthermore, ifverified = True
is also set, while the left endpoints of the intervals are lower bounds for the maximal cusp area matrix entries, the right endpoints are meaningless and could be smaller or larger than the maximal cusp area matrix entries.Verified computation
If
verified = False
, floating-point issues can arise resulting in incorrect values. The method can be made verified by passingverified = True
:sage: M.cusp_area_matrix(verified=True) # doctest: +NUMERIC3 [ 28.0000? 7.000000000000? 7.00000000000?] [7.000000000000? 28.000000? 7.00000000000?] [ 7.00000000000? 7.00000000000? 28.00000?]
- Parameters:
verified – Use verified computation.
bits_prec – Precision used for computation. Increase if computation did not succeed or a more precise result is desired.
method – Switches to older algorithms giving lower bounds when
trigDependentTryCanonize
andtrigDependent
.
- Returns:
Maximal cusp area matrix (default) or lower bounds (if
method
switches to older algorithm).
- cusp_areas(policy: str = 'unbiased', method: str = 'maximal', verified: bool = False, bits_prec: int | None = None, first_cusps: List[int] = [])
Returns a list of areas, one for each cusp. The cusp neighborhoods defined by these areas are embedded and disjoint. Furthermore, these neighborhoods are maximal in that they fail to be embedded or disjoint if any cusp neighborhood is enlarged (unless
method
is set to a value different from the default).There are different policies how these cusp neighborhoods are found.
The default
policy
isunbiased
. This means that the cusp neighborhoods are blown up simultaneously and a cusp neighborhood stops growing when it touches any cusp neighborhood including itself:>>> M = Manifold("s776") >>> M.cusp_areas() [2.64575131106459, 2.64575131106459, 2.64575131106459]
Alternatively,
policy='greedy'
can be specified. This 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, and so on:>>> M.cusp_areas(policy='greedy') [5.29150262212918, 1.32287565553230, 1.32287565553229]
Use
first_cusps
to specify the order in which the cusp neighborhoods are blown up:>>> M.cusp_areas(policy='greedy', first_cusps=[1,0,2]) [1.32287565553230, 5.29150262212918, 1.32287565553229]
An incomplete list can be given to
first_cusps
. In this case, the list is automatically completed by appending the remaining cusps in order. Thus, the above call is equivalent to:>>> M.cusp_areas(policy='greedy', first_cusps=[1]) [1.32287565553230, 5.29150262212918, 1.32287565553229]
Under the hood, this method is using
cusp_area_matrix()
.Verified computation
If
verified = False
, floating-point issues can arise resulting in incorrect values. The method can be made verified by passingverified = True
:sage: M=Manifold("s776") sage: M.cusp_areas(verified=True) # doctest: +NUMERIC9 [2.64575131107?, 2.64575131107?, 2.64575131107?]
- Parameters:
verified – Use verified computation.
bits_prec – Precision used for computation. Increase if computation did not succeed or a more precise result is desired.
method – Passed to
cusp_area_matrix()
. If set to a value different from the defaultmaximal
, the cusp neighborhoods stop growing when the corresponding value in the computed cusp area matrix is exceeded. At this point, the cusp neighborhood might not necessarily touch any other cusp neighborhood since we do not use the maximal cusp area matrix.policy – Specifies process of choosing cusp neighborhoods. Either
unbiased
orgreedy
, see above.first_cusps – Preference order of cusps. Only relevant if
policy='greedy'
, see above.
- Returns:
Areas of maximal embedded and disjoint cusp neighborhoods (default). Or areas of some embedded and disjoint cusp neighborhoods (if
method
switches to older algorithm).
- 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() CuspNeighborhood
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: str = 'unbiased', method: str = 'maximal', verified: bool = False, bits_prec: int | None = None, first_cusps: List[int] = [])
Returns a list of the (complex) Euclidean translations corresponding to the meridian and longitude of each cusp.
That is, the method uses
cusp_areas()
to find (maximal) embedded and disjoint cusp neighborhoods. It then uses the boundaries of these cusp neighborhoods to measure the meridian and longitude of each cusp. The result is a pair for each cusp. The first entry of the pair corresponds to the meridian and is complex. The second entry corresponds to the longitude and 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)]
It takes the same arguments as
cusp_areas()
:>>> M.cusp_translations(policy = 'greedy') [(0.70710678118654752440084436210 + 1.8708286933869706927918743662*I, 2.8284271247461900976033774484), (0.35355339059327376220042218105 + 0.93541434669348534639593718308*I, 1.4142135623730950488016887242), (0.35355339059327376220042218105 + 0.93541434669348534639593718308*I, 1.4142135623730950488016887242)]
Verified computations
If
verified = False
, floating-point issues can arise resulting in incorrect values. The method can be made verified by passingverified = True
:sage: M.cusp_translations(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?)]
Note that the first element of each pair is a SageMath
ComplexIntervalField
and the second element aRealIntervalField
.
- dehn_fill(filling_data, which_cusp=None) 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=1e-08, displacement=(0.0, 0.0, 0.0), centroid_at_origin: bool = True, maximize_injectivity_radius: bool = True, include_words: bool = False) DirichletDomain
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: int | None = None, verbose: bool = False) Manifold
Drills the geodesic corresponding to the given word in the unsimplified fundamental group. Here is an example:
>>> M = Manifold("m004") >>> M.length_spectrum_alt(max_len=1.2) [Length Core curve Word 1.08707014499574 + 1.72276844987009*I - bC, 1.08707014499574 - 1.72276844987009*I - a] >>> 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 resulting manifold corresponds to the drilled geodesic. The longitude and meridian for that cusp are chosen such that
(1,0)
-filling the last cusp results in the given (undrilled) manifold:>>> N.dehn_fill((1,0),-1) >>> M.is_isometric_to(N) True >>> N.cusp_info(1)['core_length'] 1.08707014499574 - 1.72276844987009*I
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 given manifold.
If the given geodesic coincides with a core curve of a filled cusp, the cusp is unfilled instead:
>>> 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
In this case, the peripheral information is also updated such that the above remark about
(1,0)
-filling applies again:>>> N.dehn_fill((1,0), -1) >>> N.volume() 1.73712388065
That is, the longitude and meridian of the unfilled cusps are reinstalled and the cusps reindexed so that the unfilled cusp becomes the last cusp.
Here is another example where we drill the core geodesic:
>>> 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]
While the result of drilling a geodesic is a triangulation and thus combinatorial in nature, some intermediate computations (for example, to compute the intersections of the geodesic with the faces of the tetrahedra) are numerical. Sometimes, it is necessary to increase the precision with
bits_prec
to make the method succeed and produce the correct result.Verified computation
If
verified = False
, floating-point issues can arise resulting in drilling the wrong loop. The method can be made verified by passingverified = True
:sage: M = Manifold("m004(2,3)") sage: M.drill_word('caa', verified = True, bits_prec = 100) m004_drilled(2,3)(0,0)
That is, if the precision is insufficient to prove the result is correct, the algorithm fails with an exception (most likely
InsufficientPrecisionError
).- Parameters:
word – The word in the unsimplified fundamental group specifying the geodesic to be drilled.
bits_prec – The precision used in the intermediate computation. Increase if the computation failed.
verified – Use verified computation.
verbose – Print intermediate results and statistics.
- Returns:
Manifold obtained by drilling geodesic.
(1,0)
-filling the last cusp gives the given (undrilled) manifold.
- drill_words(words: Sequence[str], verified: bool = False, bits_prec: int | None = None, verbose: bool = False) Manifold
A generalization of
drill_word
to drill several geodesics simultaneously. It takes a list of words in the unsimplified fundamental group.Here is an example where we drill two geodesics. One of the geodesics is the core curve corresponding to the third cusp. The other geodesic is not a core curve:
>>> 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)
Let n be the number of geodesics that were drilled. Then the last n cusps correspond to the drilled geodesics and appear in the same order than the geodesics were given as words. Note that in the above example, we expect six cusps since we started with four cusps and drilled two geodesics. However, we only obtain five cusps 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 given (undrilled) manifold by
(1,0)
-filling the last n cusps.>>> N.dehn_fill((1,0), -2) >>> N.dehn_fill((1,0), -1) >>> 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]
- Parameters:
word – The words in the unsimplified fundamental group specifying the geodesics to be drilled.
bits_prec – The precision used in the intermediate computation. Increase if the computation failed.
verified – Use verified computation.
verbose – Print intermediate results and statistics.
- Returns:
Manifold obtained by drilling geodesics.
(1,0)
-filling the last n cusps gives the given (undrilled) manifold where n is the number of given words.
- 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}
- exterior_to_link(verbose: bool = False, check_input: bool = True, check_answer: bool = True, careful_perturbation: bool = True, simplify_link: bool = True, pachner_search_tries: int = 10, seed: int | None = None) Link
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:
If no cusp is filled, then they are the meridians of the current peripheral curves.
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
: WhenTrue
, prints progress updates as the algorithm goes along.check_input
: WhenTrue
(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 passcheck_input=False
.check_answer
: WhenTrue
(the default), take the exterior of the final link diagram and useManifold.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. Whencareful_perturbation=True
(the default), computations are performed to ensure this rounding does not change the isotopy class of the link.simplify_link
: WhenTrue
(the default), usesLink.simplify('global')
to minimize the size of the final diagram; otherwise, it just doesbasic
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
andcareful_perturbation
isTrue
, 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') snappy.Manifold
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: bool = True, fillings_may_affect_generators: bool = True, minimize_number_of_generators: bool = True, try_hard_to_shorten_relators: bool = True) HolonomyGroup
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')
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
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 * z'_{0000,0} * z_{0000,1} * {z''_{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), azp
the cross ratio \(z'\) at the edge (0,2) and azpp
the cross ratio \(z''\) at the edge (1,2). The entire symbolz_xxxx_y
then denotes the cross ratio belonging to the subsimplex at integral pointxxxx
(always0000
forN = 2
) of the simplexy
.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() bool
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() AbelianGroup
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) []
- inside_view(cohomology_class=None, geodesics: Sequence[str] = [])
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, \partial M; \mathbb{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 andn
- 1 to pick thei
-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.
Geodesics can be specified as words in the unsimplified fundamental group:
>>> M = Manifold("m004") >>> M.inside_view(geodesics=['a', 'bC'])
- 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(other: Manifold | ManifoldHP, return_isometries: bool = False) bool | List[Isometry]
Returns
True
if M and N are isometric,False
if they not. ARuntimeError
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 answerFalse
may not be as there could be numerical errors resulting in finding an incorrect canonical triangulation.
- is_orientable() bool
Return whether the underlying 3-manifold is orientable.
>>> M = Triangulation('x124') >>> M.is_orientable() False
- is_two_bridge() bool
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 answerFalse
, as there could be numerical errors resulting in finding an incorrect canonical triangulation.
- isometry_signature(of_link=False, ignore_orientation=True, verified=False, interval_bits_precs=[53, 212], exact_bits_prec_and_degrees=[(212, 10), (1000, 20), (2000, 20)], verbose=False) str
Returns the “isometry signature”, a complete invariant of the hyperbolic 3-manifold obtained by applying the Dehn-fillings. The isometry signature is always a (decorated) isomorphism signature, see
triangulation_isosig()
, and was introduced in Goerner ‘16.Depending on
ignore_orientation
, it is a complete invariant of either the oriented (if orientable) or unoriented hyperbolic 3-manifold. Ifof_link = True
is specified, the signature is decorated by the unoriented peripheral curves (aka meridian and longitude, up to homotopy). If the 3-manifold arises as a link complement, the decorated isometry signature obtained withof_link = True
is a complete invariant of the link.The isometry signature is computed differently based on whether there is at least one unfilled cusp.
Cusped manifolds
If there is at least one unfilled cusped, we are in the cusped case.
Here is an example of two links having isometric (hyperbolic) complements:
>>> M = Manifold("L5a1") >>> N = Manifold("L7n2") >>> M.isometry_signature() 'eLPkbdcddhgggb' >>> N.isometry_signature() 'eLPkbdcddhgggb'
The complements do have opposite handedness though:
>>> M.isometry_signature(ignore_orientation=False) 'eLPkbdcddxvvcv' >>> N.isometry_signature(ignore_orientation=False) 'eLPkbdcddhgggb'
We can show that the two links are distinct:
>>> M.isometry_signature(of_link = True) 'eLPkbdcddhgggb_baCbbaCb' >>> N.isometry_signature(of_link = True) 'eLPkbdcddhgggb_bBcBbaCb'
If we Dehn-fill some cusps, the method uses the filled triangulation. Here, we Dehn-fill the Whitehead link to get the figure-eight knot:
>>> M.dehn_fill((1,1), 0) >>> M.isometry_signature(of_link = True) 'cPcbbbiht_bacb' >>> Manifold("4_1").isometry_signature(of_link = True) 'cPcbbbiht_bacb'
In general, the isometry signature is the isomorphism signature (see
triangulation_isosig()
) of thecanonical_retriangulation()
of thefilled_triangulation()
:>>> T = M.filled_triangulation().canonical_retriangulation() >>> T.triangulation_isosig(ignore_cusp_ordering = True, ... ignore_curve_orientations = True) 'cPcbbbiht_bacb'
Closed manifolds
If all cusps are filled, we are in the closed case. In this case, the isometry signature gives the resulting closed hyperbolic 3-manifold as canonical surgery on a hyperbolic 1-cusped manifold (which is encoded by its isometry signature). Only orientable manifolds are supported in the closed case.
>>> M = Manifold("v2000(1,3)") >>> M.isometry_signature() 'fLLQcacdedenbxxrr(-7,12)'
The following code illustrates how the isometry signature is computed:
>>> M.length_spectrum_alt(count=2) [Length Core curve Word 0.06491027903143 - 2.63765810995071*I - d, 0.49405010583448 + 2.38451103485706*I - a] >>> K = M.drill_word('d').filled_triangulation().canonical_retriangulation() >>> K.dehn_fill((1,0), 0) >>> K.triangulation_isosig(ignore_cusp_ordering=True, ignore_curves=True) 'fLLQcacdedenbxxrr(-7,12)'
Note that there is clearly a unique shortest geodesic in this example. In general, the method first considers a canonical set of geodesics. For each such geodesic, it computes a candidate signature as above. It then picks a canonical signature among the candidates. Further details can be found in an upcoming paper.
Verified computations
While the isometry signature is purely combinatorial, some intermediate computations are numerical. Thus, if
verified = False
, floating-point issues can arise.The method can be made verified by passing
verified = True
:sage: M=Manifold("m007(4,1)") sage: M.isometry_signature(verified=True) 'eLPkbcdddhggsj(3,1)'
This method always needs to compute at least one canonical retriangulation. It can take the same arguments as
canonical_retriangulation()
and passes them tocanonical_retriangulation()
when computing the verified canonical retriangulation. If the manifold is closed, interval arithmetic is used when finding and drilling the short geodesics.- Parameters:
of_link – Also encode the unoriented peripheral curves. Note that it is not necessary for the manifold to be a link complement to invoke this flag. Only relevant in the cusped case.
ignore_orientation – Do not encode the orientation of the 3-manifold.
verified – Use verified computation.
interval_bits_precs – Passed to
canonical_retriangulation()
and (in the closed case) also used when callinglength_spectrum_alt_gen()
anddrill_word()
to find and drill the short geodesics.exact_bits_prec_and_degrees – Passed to
canonical_retriangulation()
.verbose – Print information about finding and drilling the short geodesics. Also passed to
canonical_retriangulation()
.
- isomorphisms_to(other: Triangulation | TriangulationHP) List[Isometry]
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: bool = True, grouped: bool = True, include_words: bool = 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
- length_spectrum_alt(count: int | None = None, max_len: Any | None = None, bits_prec: int | None = None, verified: bool = False) List[LengthSpectrumGeodesicInfo]
Returns a list of geodesics. How far this list goes can be specified by either a cut-off length or a count. The method only supports orientable manifolds. It is a convenience method for
length_spectrum_alt_gen()
. We refer the reader tolength_spectrum_alt_gen()
for further details not covered here.Cut-off length
Here is an example where a cut-off length for the geodesics is specified:
>>> M = Manifold("m202(3,4)(3,4)") >>> M.length_spectrum_alt(max_len = 0.5) [Length Core curve Word 0.14820741547094 - 1.76955170166922*I Cusp 1 bcDc, 0.14820741547097 - 1.76955170166923*I Cusp 0 aabcDabcB]
It also supports verified computations:
sage: M.length_spectrum_alt(max_len = 0.5, verified = True, bits_prec = 100) # doctest: +SKIP [Length Core curve Word 0.148207415470948? - 1.76955170166924? *I Cusp 0 aabcDabcB, 0.14820741547094... - 1.76955170166923...*I Cusp 1 bcDc]
If
verified=True
, the returned list is guaranteed to include all geodesics up to the given cut-off length and might include additional geodesics.Count
Here is an example where a count is specified:
>>> M = Manifold("m202(3,4)(3,4)") >>> M.length_spectrum_alt(count = 3) [Length Core curve Word 0.14820741547094 - 1.76955170166922*I Cusp 1 bcDc, 0.14820741547097 - 1.76955170166923*I Cusp 0 aabcDabcB, 0.79356651781096 + 2.65902431489655*I - aB, 0.79356651781096 + 2.65902431489655*I - b]
Note that the number of geodesics listed might be larger than the given count. In particular, this happens when the same (real) length appears multiple times. If
verified=True
, the returned list is guaranteed to include thecount
shortest geodesics and might include additional geodesics.Verified systole
Even though, the first reported geodesic might not be the shortest, we obtain an interval containing the systole as follows, also see
length_spectrum_alt_gen()
:sage: M = Manifold("m004") sage: M.length_spectrum_alt(count=1, verified=True, bits_prec=100)[0].length.real() # doctest: +NUMERIC21 1.0870701449957390997853?
- Parameters:
count – Number of shortest geodesics to list. The actual result might contain additional geodesics. Exactly one of
count
andmax_len
has to be specified.max_len – Cut-off length for geodesics. The actual result includes all geodesics up to the given length and might include additional geodesics. Exactly one of
count
andmax_len
has to be specified.bits_prec – Precision used for the computation. Increase if computation did not succeed.
verified – Use verified computation.
- Returns:
A list of geodesics such that the (lower bound of) the real length is non-decreasing.
- length_spectrum_alt_gen(bits_prec: int | None = None, verified: bool = False) Sequence[LengthSpectrumGeodesicInfo]
Returns a generator for the geodesics sorted by real length. The method only supports orientable manifolds.
Here is an example:
>>> M = Manifold("m202(3,4)(0,0)") >>> spec = M.length_spectrum_alt_gen() >>> next(spec) Length Core curve Word 0.14742465268512 - 1.78287093565202*I Cusp 0 aabcDabcB >>> next(spec) 0.81161414965958 + 2.72911699294426*I - b >>> next(spec) 0.84163270359334 + 2.61245944742151*I - aB
Note that the shortest geodesic in the above example happens to be the core curve of the filled cusp (Cusp 0).
Access just the length or word:
>>> g = next(spec) >>> g.length 0.93461379591349 + 2.70060614107722*I >>> g.word 'a'
The word is given with respect to the unsimplified fundamental group:
>>> G = M.fundamental_group(simplify_presentation=False) >>> G.complex_length('a') 0.93461379591349 + 2.70060614107722*I
The method also supports higher precision:
>>> M = Manifold("m003(-3,1)") >>> spec = M.length_spectrum_alt_gen(bits_prec=100) >>> next(spec).length 0.58460368501798696932015666264 + 2.4953704555604684110903962008*I
Performance
This method uses a different algorithm than
length_spectrum
. In particular, it does not compute the Dirichlet domain. It also allows for verified computations. It is implemented in python and thus typically slower thanlength_spectrum
. But there are also some cases where it is significantly faster. In particular, this applies to spun triangulations such asm004(21,10)
.Here is example where we can help the algorithm by guessing and drilling and filling a short geodesic:
>>> M = Manifold("o9_00639")
Over an hour to compute:
>>> M.high_precision().length_spectrum(0.1) mult length topology parity 1 0.00150226276052 - 2.39996262244127*I circle +
A couple of minutes to compute:
>>> spec = M.length_spectrum_alt_gen(bits_prec = 150) >>> next(spec) Length Word Core curve 0.00150226276052 - 2.39996262244127*I a -
After drilling and filling, less than a second to compute:
>>> N = M.drill_word('a') >>> N.dehn_fill((1,0),-1) # N is now isometric to o9_00639 but as a surgery m125(0,0)(34,55) >>> spec = N.length_spectrum_alt_gen() >>> next(spec) Length Core curve Word 0.00150226276073 - 2.39996262244128*I Cusp 1 cDcDDcDcDDcDDcDcDDcDcDDcDDcDcDDcDD >>> next(spec).length.real() 0.96218768626877
Verified computations
The method also supports verified computations:
sage: M = Manifold("m019") sage: spec = M.length_spectrum_alt_gen(verified=True, bits_prec=100) sage: next(spec) Length Core curve Word 0.43153441294719... + 2.35105908147863...*I - a sage: next(spec) 0.88944299721255... - 2.94185904702273...*I - bD
If
verified = True
is passed, the algorithm guarantees that the lower bound of the real length is (non-strictly) increasing. In particular, we know that we have found all geodesics less than the following length:sage: next(spec).length.real().lower() # doctest: +NUMERIC12 0.94135129037387168886341739832
To illustrate some pitfalls, here is an example of a potential a result of the method:
Real length interval
Word
[1.0, 2.0]
a
[1.2, 1.3]
b
[1.7, 1.8]
c
[3.0, 4.0]
d
Note that we cannot say whether geodesic
a
is actually the first, second or third shortest geodesic or tied withb
orc
. Increasing precision can change (representative words and) the order in which the geodesics are emitted.We can say that together
a
,b
andc
are the three shortest geodesics. Furthermore, we can also say that the systole of the manifold is in[1.0, 2.0]
even thougha
itself might not be the shortest geodesic. The latter is true in general:Verified systole
It is not necessarily true that the first geodesic returned by the method is the shortest geodesic. Despite this, the interval for the real length of the first geodesic always contains the systole of the manifold:
sage: M = Manifold("m004") sage: spec = M.length_spectrum_alt_gen(verified=True) sage: g = next(spec) # g might or might not be shortest geodesic sage: systole = g.length.real() # But interval is large enough to contain systole sage: systole # doctest: +NUMERIC6 1.08707015?
- Parameters:
bits_prec – Precision used for the computation. Increase if computation did not succeed.
verified – Use verified computation.
- Returns:
A generator to enumerate the geodesics such that the (lower bound of the) real length is non-decreasing.
- link()
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()
andplink
.
- name() str
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 Dunfield and Garoufalidis ‘12 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') int
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() int
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'
- plink()
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 equipped 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)
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)
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() None
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: str) None
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)
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: str = 'unbiased', method: str = 'maximal', verified: bool = False, bits_prec: int | None = None, first_cusps: List[int] = [])
Returns a list of short slopes (for Dehn-fillings) for each cusp.
That is, the method uses
cusp_areas()
to find (maximal) embedded and disjoint cusp neighborhoods. It uses the boundaries of these cusp neighborhoods to measure the length of a peripheral curve. For each cusp, it determines all simple peripheral curves shorter than the givenlength
(which defaults to 6). The result is a list of the corresponding slopes for each cusp:>>> M = Manifold("otet20_00022") >>> M.short_slopes() [[(1, 0), (-1, 1), (0, 1)], [(1, 0)]]
It takes the same arguments as
cusp_areas()
:>>> M.short_slopes(policy = 'greedy') [[(1, 0)], [(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)]]
Verified computation
If
verified = False
, floating-point issues can arise resulting in incorrect values. The method can be made verified by passingverified = True
:sage: M = Manifold("4_1") sage: M.short_slopes(verified = True) [[(1, 0), (-4, 1), (-3, 1), (-2, 1), (-1, 1), (0, 1), (1, 1), (2, 1), (3, 1), (4, 1)]]
If
verified = True
, the result is guaranteed to contain all short slopes and might contain additional slopes (with lengths slightly longer than the givenlength
but this could not be proven using the interval estimates).The given
length
is cast to a SageMathRealIntervalField
of the given precision ifverified = True
: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, 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, 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: bool = False) SymmetryGroup
Returns the symmetry group of the Manifold. If the flag “of_link” is set, then it only returns symmetries that preserves the meridians.
- symplectic_basis(verify: bool = False)
Extend the Neumann-Zagier matrix to one which is symplectic (up to factors of 2) using oscillating curves, see Mathews and Purcell ‘22. Only accepts triangulations with 1 cusp.
>>> M = Manifold("4_1") >>> M.symplectic_basis() [-1 0 -1 -1] [ 2 0 -2 0] [-2 -1 -2 -1] [ 0 -1 -2 -1]
- Parameters:
verify – Explicitly test if the resulting matrix is symplectic.
- 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: bool = True, ignore_cusp_ordering: bool = False, ignore_curves: bool = False, ignore_curve_orientations: bool = False, ignore_filling_orientations: bool = False, ignore_orientation: bool = True) str
Returns the “(decorated) isomorphism signature”, a compact text representation of the triangulation:
>>> T = Triangulation('m004') >>> T.triangulation_isosig() 'cPcbbbiht_BaCB'
This string can be used later to recreate an isomorphic triangulation:
>>> U = Triangulation('cPcbbbiht_BaCB') >>> T == U True
The isomorphism signature is also used to compute the
isometry_signature
. It comes in two flavors controlled by thedecorated
flag.Undecorated isomorphism signature
The undecorated isomorphism signature is a complete invariant of the (oriented) triangulation up to combinatorial isomorphism:
>>> T = Triangulation('m015') >>> T.triangulation_isosig(decorated=False) 'dLQbcccdero'
It was introduced in Burton ‘11. It canonizes and generalizes the ealier dehydration string by Callahan, Hildebrand and Weeks ‘99. The undecorated isomorphism signature can also be given to Regina’s
Triangulation3.fromIsoSig
.By default, the orientation (if orientable) is ignored. More precisely, it computes the string for both orientations (if orientable) and uses the lexicographically smaller string:
>>> T = Triangulation('m015') >>> T.triangulation_isosig(decorated=False) 'dLQbcccdero' >>> T.reverse_orientation() >>> T.triangulation_isosig(decorated=False) 'dLQbcccdero'
When specifying
ignore_orientation = False
, the result encodes the orientation (if orientable). Now the result is different if we change the orientation of a chiral triangulation:>>> T = Triangulation('m015') >>> T.triangulation_isosig(decorated=False, ignore_orientation=False) 'dLQbcccdero' >>> T.reverse_orientation() >>> T.triangulation_isosig(decorated=False, ignore_orientation=False) 'dLQbccceekg'
Decorated isomorphism signature (default)
SnapPy can decorate the isomorphism signature to include the following peripheral information in a canonical way (that is invariant under the action by combinatorial isomorphisms of the triangulation):
Indexing of the cusps (that is, ideal vertices).
Included by default. Can be suppressed with
ignore_cusp_ordering = True
.
Peripheral curves (aka meridian and longitude, up to homotopy).
Included by default. Can be suppressed with
ignore_curves = True
.By default, the decoration encodes the oriented peripheral curves. By specifying
ignore_curve_orientations = True
, it encodes the unoriented peripheral curves instead.
Dehn-fillings (if present).
By default, the decoration encodes the oriented Dehn-fillings. That is, we also encodes the orientation of the peripheral curve that is used for the Dehn-filling (this explanation only works if the coefficients are integral). By specifying
ignore_filling_orientations = True
, the decoration encodes the unoriented Dehn-fillings. That is, it normalizes the Dehn-filling coefficients by picking a canonical pair among \((m,l)\) and \((-m,-l)\).
Details of the encoding are explained in the SnapPy source code.
Example
Let us consider the links \(9^2_{34}\) and
L9a21
. Note that we usecanonical_retriangulation
to make the following examples say something intrinsic about the hyperbolic manifold:>>> from snappy import Manifold >>> M = Manifold('9^2_34').canonical_retriangulation() >>> N = Manifold('L9a21').canonical_retriangulation()
The decorated isosig recovers the entire peripheral information faithfully (including orientation, see below):
>>> M.triangulation_isosig() 'oLLvzQLLQQccdhifihnlmkmlnnpvuvbvouggbggoo_baBabbbBbC' >>> K = Triangulation('oLLvzQLLQQccdhifihnlmkmlnnpvuvbvouggbggoo_baBabbbBbC') >>> M.isomorphisms_to(K) [0 -> 0 1 -> 1 [1 0] [1 0] [0 1] [0 1] Extends to link]
The two links have isometric complements:
>>> M.triangulation_isosig(decorated=False) 'oLLvzQLLQQccdhifihnlmkmlnnpvuvbvouggbggoo' >>> N.triangulation_isosig(decorated=False) 'oLLvzQLLQQccdhifihnlmkmlnnpvuvbvouggbggoo'
However, the complements have different handedness:
>>> M.triangulation_isosig(decorated=False,ignore_orientation=False) 'oLLzLPwzQQccdeghjiiklmnmnnuvuvvavovvffffo' >>> N.triangulation_isosig(decorated=False,ignore_orientation=False) 'oLLvzQLLQQccdhifihnlmkmlnnpvuvbvouggbggoo'
Also, the cusps/components of the link are indexed differently:
>>> M.triangulation_isosig(ignore_curves=True) 'oLLvzQLLQQccdhifihnlmkmlnnpvuvbvouggbggoo_ba' >>> N.triangulation_isosig(ignore_curves=True) 'oLLvzQLLQQccdhifihnlmkmlnnpvuvbvouggbggoo_ab'
Ignoring the indexing, we also see that the oriented merdians and longitudes do not match:
>>> M.triangulation_isosig(ignore_cusp_ordering=True) 'oLLvzQLLQQccdhifihnlmkmlnnpvuvbvouggbggoo_bBbCBabb' >>> N.triangulation_isosig(ignore_cusp_ordering=True) 'oLLvzQLLQQccdhifihnlmkmlnnpvuvbvouggbggoo_BbbCbabb'
However, they are the same links (ignoring indexing and orientation):
>>> M.triangulation_isosig(ignore_cusp_ordering=True, ignore_curve_orientations=True) 'oLLvzQLLQQccdhifihnlmkmlnnpvuvbvouggbggoo_bBBcbabb' >>> N.triangulation_isosig(ignore_cusp_ordering=True, ignore_curve_orientations=True) 'oLLvzQLLQQccdhifihnlmkmlnnpvuvbvouggbggoo_bBBcbabb'
Let us create two surgery presentations from the links (note that we fill after
canonical_retriangulation
since it rejects Dehn-fillings):>>> M.dehn_fill((4, 5),0) >>> N.dehn_fill((4,-5),1)
They are equivalent surgery presentations (of the same manifold):
>>> M.triangulation_isosig( ... ignore_cusp_ordering=True, ... ignore_curves=True, ... ignore_filling_orientations=True) 'oLLvzQLLQQccdhifihnlmkmlnnpvuvbvouggbggoo(0,0)(1,5)' >>> N.triangulation_isosig( ... ignore_cusp_ordering=True, ... ignore_curves=True, ... ignore_filling_orientations=True) 'oLLvzQLLQQccdhifihnlmkmlnnpvuvbvouggbggoo(0,0)(1,5)'
Orientation
Note that
ignore_orientation=True
only applies to the undecorated part of the isomorphism signature. The decoration can still capture the the orientation. More, precisely, the result oftriangulation_isosig()
depends on the orientation (if the triangulation is orientable and chiral) if any of the following is true:ignore_orientation = False
.decorated = True
andignore_curves = False
andignore_filling_orientations = False
.
In these cases, re-constructing a triangulation from the isomorphism signature yields a triangulation with the same handedness.
- Parameters:
decorated – Include peripheral information such as indexing of the cusps, (oriented or unoriented) peripheral curves and (oriented or unoriented) Dehn-fillings.
ignore_cusp_ordering – Do not encode the indexing of the cusps. Only relevant if
decorated = True
.ignore_curves – Do not encode the peripheral curves. Only relevant if
decorated = True
. This is new in SnapPy version 3.2. Ifignore_curves = True
, the result of this method cannot be given to prior versions.ignore_curve_orientations – Do not encode the orientations of the peripheral curves. Only relevant if
decorated = True
andignore_curves = False
.ignore_filling_orientations – Do not encode the orientations of the Dehn-fillings. Only relevant if
decorated = True
.ignore_orientation – Do not encode the orientation of the triangulation in the undecorated part of the triangulation isosig. See above section about orientation.
- classmethod use_field_conversion(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: R = M.volume().parent() sage: R.name().startswith('RealField') True sage: from snappy.number import SnapPyNumbers sage: Manifold.use_field_conversion('snappy') sage: M = Manifold('m004') sage: R = M.volume().parent() sage: isinstance(R, SnapPyNumbers) True sage: Manifold.use_field_conversion('sage') sage: M = Manifold('m004') sage: R = M.volume().parent() sage: R.name().startswith('RealField') True
- 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’sComplexIntervalField
) 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 callscheck_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 aManifold
.>>> 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