Classes

PtolemyVariety

class snappy.ptolemy.ptolemyVariety.PtolemyVariety(manifold, N, obstruction_class, simplify, eliminate_fixed_ptolemys)

Holds a reduced Ptolemy variety.

=== Examples ===

To generate such a variety, call:

>>> from snappy import Manifold
>>> p = Manifold("4_1").ptolemy_variety(2, obstruction_class = 1)

Show the equations and variables:

>>> for e in p.equations: print(e)
- 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
>>> p.variables
['c_0011_0', 'c_0101_0']

Show as an ideal (sage object):

>>> p.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
(skip doctest because example only works in sage and not plain python)

Produce magma input:

>>> s = p.to_magma()
>>> print(s.split('ring and ideal')[1].strip())    
R<c_0011_0, c_0101_0> := PolynomialRing(RationalField(), 2, "grevlex");
MyIdeal := ideal<R |
          - c_0011_0 * c_0101_0 + c_0011_0^2 + c_0101_0^2,
    ...

Call p.compute_solutions() to automatically compute solutions!

Show canonical representatives:

(The Ptolemy coordinates c_0110_0 and c_0101_0 are identified, this information is needed to recover all Ptolemy coordinates from the solutions of a simplified Ptolemy variety. The information is also packaged into a python section by py_eval_variable_dict().)

>>> p.canonical_representative["c_0110_0"]
(1, 0, 'c_0101_0')
compute_decomposition(engine=None, memory_limit=750000000, directory=None, verbose=False, template_path='magma/default.magma_template')

Starts an engine such as magma to compute the radical decomposition of the Ptolemy variety.

If started in sage, uses sage, otherwise, uses magma.

=== Arguments ===

engine — engine to use, currently, only support magma and sage memory_limit — maximal allowed memory in bytes directory — location for input and output files, temporary directory used if not specified verbose — print extra information

compute_solutions(engine=None, numerical=False, template_path='magma/default.magma_template', memory_limit=750000000, directory=None, verbose=False)

Starts an engine such as magma to compute the radical decomposition of the ideal and computes exact solutions.

If started in sage, uses sage, otherwise, uses magma.

=== Arguments ===

engine — engine to use, currently, only support magma and sage numerical — get numerical solutions from magma instead of exact ones memory_limit — maximal allowed memory in bytes directory — location for input and output files, temporary directory used if not specified verbose — print extra information

degree_to_shapes()

In general, there can be d different solutions to the (reduced) Ptolemy variety for each solution to the gluing equations (with peripheral equations). This method computes d which is also the number of elements in H^1(Mhat; Z/N) where Mhat is the non-manifold cell comples obtained by gluing together the tetrahedra as non-ideal tetrahedra.

For example, the Ptolemy variety for m009 has 4 points but there are only 2 distinct boundary-unipotent PSL(2,C) representations. Thus the following call returns 2 = 4 / 2

>>> from snappy import Manifold
>>> Manifold("m009").ptolemy_variety(2,1).degree_to_shapes()
2
>>> Manifold("m010").ptolemy_variety(2,1).degree_to_shapes()
2
>>> Manifold("m011").ptolemy_variety(2,1).degree_to_shapes()
1
>>> Manifold("m009").ptolemy_variety(3,1).degree_to_shapes()
1
>>> Manifold("m010").ptolemy_variety(3,1).degree_to_shapes()
3
>>> Manifold("m011").ptolemy_variety(3,1).degree_to_shapes()
1
filename_base()

Preferred filename base for writing out this Ptolemy variety

>>> from snappy import *
>>> M = Manifold('4_1')
>>> p = M.ptolemy_variety(2, obstruction_class = 1)
>>> p.filename_base()
'4_1__sl2_c1'
>>> p = M.ptolemy_variety(2)
>>> p.filename_base()
'4_1__sl2_c0'
py_eval_section()

Returns a string that can be evaluated in python and contains extra information needed to recover solutions from a simplified Ptolemy variety.

>>> from snappy import Manifold, pari
>>> M = Manifold('4_1')
>>> p = M.ptolemy_variety(2, obstruction_class = 1)

Get extra information

>>> eval_section = p.py_eval_section()
>>> print(eval_section)    
{'variable_dict' :
     (lambda d: {
    ...

Turn it into a python object by evaluation.

>>> obj = eval(eval_section)

Access the function that expands a solution to the simplified Ptolemy variety to a full solution.

>>> variable_dict = obj['variable_dict']

Setup a solution and expand it to a full solution, ‘1’ must map to 1

>>> simplified_solution = {'c_0101_0' : pari('0.5 - 0.866025403784439*I'), '1' : pari(1), 'c_0011_0' : pari(1)}
>>> full_solution = variable_dict(simplified_solution)

Full solution is a dictionary with a key for every Ptolemy coordinate

>>> full_solution['c_1010_1']
1
>>> for tet in range(2):
...     for i in utilities.quadruples_with_fixed_sum_iterator(2, skipVertices = True):
...         c = "c_%d%d%d%d" % i + "_%d" % tet
...         assert c in full_solution
to_magma(template_path='magma/default.magma_template')

Returns a string with the ideal that can be used as input for magma.

The advanced user can provide a template string to write own magma code to process the equations.

>>> from snappy import *
>>> p = Manifold("4_1").ptolemy_variety(2, obstruction_class = 1)

Magma input to compute radical Decomposition >>> s = p.to_magma() >>> print(s.split(‘ring and ideal’)[1].strip()) #doctest: +ELLIPSIS +NORMALIZE_WHITESPACE R<c_0011_0, c_0101_0> := PolynomialRing(RationalField(), 2, “grevlex”); MyIdeal := ideal<R | - c_0011_0 * c_0101_0 + c_0011_0^2 + c_0101_0^2, … >>> “RadicalDecomposition” in p.to_magma() True

to_magma_file(filename, template_path='magma/default.magma_template')
>>> import os, tempfile
>>> from snappy import Manifold
>>> handle, name = tempfile.mkstemp()
>>> p = Manifold("4_1").ptolemy_variety(2, obstruction_class=1)
>>> p.to_magma_file(name)
>>> os.close(handle); os.remove(name)

PtolemyCoordinates

class snappy.ptolemy.coordinates.PtolemyCoordinates(d, is_numerical=True, py_eval_section=None, manifold_thunk=<function PtolemyCoordinates.<lambda>>, non_trivial_generalized_obstruction_class=False)

Represents a solution of a Ptolemy variety as python dictionary.

=== Examples ===

Construct solution from magma output:

>>> from snappy.ptolemy.processMagmaFile import _magma_output_for_4_1__sl3, solutions_from_magma
>>> from snappy import Manifold
>>> solutions = solutions_from_magma(_magma_output_for_4_1__sl3)
>>> solution = solutions[2]

Access a Ptolemy coordinate:

>>> solution['c_2100_0']
1
>>> solution.number_field()
x^2 + x + 1

Solution is always 0 dimensional:

>>> solution.dimension
0

Check that it is really a solution, exactly:

>>> solution.check_against_manifold()

If the solution was not created through the ptolemy module and thus not associated to a manifold, we need to explicitly specify one:

>>> myDict = {}
>>> for key, value in solution.items():
...     myDict[key] = value
>>> mysolution = PtolemyCoordinates(myDict)
>>> M = Manifold("4_1")
>>> solution.check_against_manifold(M)

Turn into (Galois conjugate) numerical solutions:

>>> old_precision = pari.set_real_precision(100) # with high precision
>>> numerical_solutions = solution.numerical()

Check that it is a solution, numerically:

>>> numerical_solutions[0].check_against_manifold(M, 1e-80)
>>> pari.set_real_precision(old_precision) # reset pari engine
100

Compute cross ratios from the Ptolemy coordinates (cross ratios according to SnapPy convention, see help(solution.cross_ratios):

>>> cross = solution.cross_ratios()
>>> cross['z_0001_0']
Mod(-x, x^2 + x + 1)

Compute volumes:

>>> volumes = cross.volume_numerical()

Check that volume is 4 times the geometric one:

>>> volume = volumes[0].abs()
>>> diff = abs(4 * M.volume() - volume)
>>> diff < 1e-9
True

Compute flattenings:

>>> flattenings = solution.flattenings_numerical()

Compute complex volumes:

>>> cvols = [flattening.complex_volume() for flattening in flattenings]
>>> volume = cvols[0].real().abs()
>>> chernSimons = cvols[0].imag()
>>> diff = abs(4 * M.volume() - volume)
>>> diff < 1e-9
True
>>> from snappy import pari
>>> normalized = chernSimons * 6 / (pari('Pi')**2)

Check that Chern Simons is zero up to 6 torsion:

>>> normalized - normalized.round() < 1e-9
True
N()

Get the N where these coordinates are for SL/PSL(N, C)-representations.

check_against_manifold(M=None, epsilon=None)

Checks that the given solution really is a solution to the Ptolemy variety of a manifold. See help(ptolemy.PtolemyCoordinates) for more example.

=== Arguments ===

  • M — manifold to check this for

  • epsilon — maximal allowed error when checking the relations, use None for exact comparison.

complex_volume_numerical(drop_negative_vols=False, with_modulo=False)

Turn into (Galois conjugate) numerical solutions and compute complex volumes. If already numerical, return the volume.

Complex volume is defined up to i*pi**2/6.

See numerical(). If drop_negative_vols = True is given as optional argument, only return complex volumes with non-negative real part.

cross_ratios()

Compute cross ratios from Ptolemy coordinates. The cross ratios are according to the SnapPy convention, so we have:

z = 1 - 1/zp, zp = 1 - 1/zpp, zpp = 1 - 1/z

where:

z   is at the edge 01 and equal to   s0 * s1 * (c_1010 * c_0101) / (c_1001 * c_0110)
zp  is at the edge 02 and equal to - s0 * s2 * (c_1001 * c_0110) / (c_1100 * c_0011)
zpp is at the edge 03 and equal to   s0 * s3 * (c_1100 * c_0011) / (c_0101 * c_1010).

Note that this is different from the convention used in Garoufalidis, Goerner, Zickert: Gluing Equations for PGL(n,C)-Representations of 3-Manifolds http://arxiv.org/abs/1207.6711

Take an exact solution:

>>> from snappy.ptolemy.processMagmaFile import _magma_output_for_4_1__sl3, solutions_from_magma
>>> solutions = solutions_from_magma(_magma_output_for_4_1__sl3)
>>> solution = solutions[2]

Turn into cross Ratios:

>>> crossRatios = solution.cross_ratios()

Get a cross ratio:

>>> crossRatios['zp_0010_0']
Mod(-x, x^2 + x + 1)

Check the relationship between cross ratios:

>>> crossRatios['z_0010_0'] == 1 - 1 / crossRatios['zp_0010_0']
True
>>> crossRatios['zp_0010_0'] == 1 - 1 / crossRatios['zpp_0010_0']
True
>>> crossRatios['zpp_0010_0'] == 1 - 1 / crossRatios['z_0010_0']
True

Get information about what one can do with cross ratios

cross_ratios_numerical()

Turn exact solutions into numerical and then compute cross ratios. See numerical() and cross_ratios().

diamond_coordinate(tet, v0, v1, v2, pt)

Returns the diamond coordinate for tetrahedron with index tet for the face with vertices v0, v1, v2 (integers between 0 and 3) and integral point pt (quadruple adding up to N-2).

See Definition 10.1: Garoufalidis, Goerner, Zickert: Gluing Equations for PGL(n,C)-Representations of 3-Manifolds http://arxiv.org/abs/1207.6711

evaluate_word(word, G=None)

Given a word in the generators of the fundamental group, compute the corresponding matrix. By default, these are the generators of the unsimplified presentation of the fundamental group. An optional SnapPy fundamental group can be given if the words are in generators of a different presentation, e.g., c.evaluate_word(word, M.fundamental_group(True)) to evaluate a word in the simplified presentation returned by M.fundamental_group(True).

For now, the matrix is returned as list of lists.

flattenings_numerical()

Turn into numerical solutions and compute flattenings, see help(snappy.ptolemy.coordinates.Flattenings) Also see numerical()

Get Ptolemy coordinates.

>>> from snappy.ptolemy.processMagmaFile import _magma_output_for_4_1__sl3, solutions_from_magma
>>> solutions = solutions_from_magma(_magma_output_for_4_1__sl3)
>>> solution = solutions[2]

Compute a numerical solution

>>> flattenings = solution.flattenings_numerical()

Get more information with help(flattenings[0])

get_manifold()

Get the manifold for which this structure represents a solution to the Ptolemy variety.

has_obstruction()

Whether the Ptolemy variety has legacy obstruction class that modifies the Ptolemy relation to

is_geometric(epsilon=1e-06)

Returns true if all shapes corresponding to this solution have positive imaginary part.

If the solutions are exact, it returns true if one of the corresponding numerical solutions is geometric.

An optional epsilon can be given. An imaginary part of a shape is considered positive if it is larger than this epsilon.

long_edge(tet, v0, v1, v2)

The matrix that labels a long edge from v0 to v1 (integers between 0 and 3) of a (doubly) truncated simplex corresponding to an ideal tetrahedron with index tet.

This matrix was labeled alpha^{v0v1v2} (v2 does not matter for non double-truncated simplex) in Figure 18 of Garoufalidis, Goerner, Zickert: Gluing Equations for PGL(n,C)-Representations of 3-Manifolds http://arxiv.org/abs/1207.6711

It is computed using equation 10.4. Note that the ratio coordinate is different from the definition in the paper (see ratio_coordinate).

The resulting matrix is given as a python list of lists.

middle_edge(tet, v0, v1, v2)

The matrix that labels a middle edge on the face v0, v1, v2 (integers between 0 and 3) of a doubly truncated simplex (or a short edge of the truncated simplex) corresponding to an ideal tetrahedron with index tet.

This matrix was labeled beta^{v0v1v2} in Figure 18 of Garoufalidis, Goerner, Zickert: Gluing Equations for PGL(n,C)-Representations of 3-Manifolds http://arxiv.org/abs/1207.6711

It is computed using equation 10.4.

The resulting matrix is given as a python list of lists.

multiply_and_simplify_terms_in_RUR()

If a Ptolemy coordinate is given as Rational Univariate Representation with numerator and denominator being a product, multiply the terms, reduce the fraction and return the result.

See multiply_and_simplify_terms of RUR.

This loses information about how the numerator and denominator are factorised.

multiply_terms_in_RUR()

If a Ptolemy coordinate is given as Rational Univariate Representation with numerator and denominator being a product, multiply the terms and return the result.

See multiply_terms of RUR.

This loses information about how the numerator and denominator are factorised.

num_tetrahedra()

The number of tetrahedra for which we have Ptolemy coordinates.

number_field()

For an exact solution, return the number_field spanned by the Ptolemy coordinates. If number_field is Q, return None.

numerical()

Turn exact solutions into numerical solutions using pari.

Take an exact solution:

>>> from snappy.ptolemy.processMagmaFile import _magma_output_for_4_1__sl3, solutions_from_magma
>>> solutions = solutions_from_magma(_magma_output_for_4_1__sl3)
>>> solution = solutions[2]

Turn into a numerical solution:

>>> old_precision = pari.set_real_precision(100) # with high precision
>>> numerical_solutions = solution.numerical()
>>> pari.set_real_precision(old_precision) # reset pari engine
100

Pick one of the Galois conjugates:

>>> numerical_solution = numerical_solutions[0]
>>> value = numerical_solution['c_1110_0']
ratio_coordinate(tet, v0, v1, pt)

Returns the ratio coordinate for tetrahedron with index tet for the edge from v0 to v1 (integers between 0 and 3) and integral point pt (quadruple adding up N-1) on the edge.

See Definition 10.2: Garoufalidis, Goerner, Zickert: Gluing Equations for PGL(n,C)-Representations of 3-Manifolds http://arxiv.org/abs/1207.6711

Note that this definition turned out to have the wrong sign. Multiply the result by -1 if v1 < v0 and N is even.

short_edge(tet, v0, v1, v2)

Returns the identity. This is because we can think of the matrices given by Ptolemy coordinates of living on truncated simplices which can be though of as doubly truncated simplices where all short edges are collapsed, hence labeled by the identity.

See equation 10.4 in Garoufalidis, Goerner, Zickert: Gluing Equations for PGL(n,C)-Representations of 3-Manifolds http://arxiv.org/abs/1207.6711

to_PUR()

If any Ptolemy coordinates are given as Rational Univariate Representation, convert them to Polynomial Univariate Representation and return the result.

See to_PUR of RUR.

This conversion might lead to very large coefficients.

volume_numerical(drop_negative_vols=False)

Turn into (Galois conjugate) numerical solutions and compute volumes. If already numerical, only return the one volume. See numerical().

If drop_negative_vols = True is given as optional argument, only return non-negative volumes.

CrossRatios

class snappy.ptolemy.coordinates.CrossRatios(d, is_numerical=True, manifold_thunk=None)

Represents assigned shape parameters/cross ratios as dictionary. The cross ratios are according to SnapPy convention, so we have:

z = 1 - 1/zp, zp = 1 - 1/zpp, zpp = 1 - 1/z

where:

z   is at the edge 01 and equal to s0 * s1 * (c_1010 * c_0101) / (c_1001 * c_0110)
zp  is at the edge 02 and equal to s0 * s2 * (c_1001 * c_0110) / (c_1100 * c_0011)
zpp is at the edge 03 and equal to s0 * s3 * (c_1100 * c_0011) / (c_0101 * c_1010).

Note that this is different from the convention used in Garoufalidis, Goerner, Zickert: Gluing Equations for PGL(n,C)-Representations of 3-Manifolds http://arxiv.org/abs/1207.6711

N()

Get the N such that these cross ratios are for SL/PSL(N,C)-representations.

check_against_manifold(M=None, epsilon=None)

Checks that the given solution really is a solution to the PGL(N,C) gluing equations of a manifold. Usage similar to check_against_manifold of PtolemyCoordinates. See help(ptolemy.PtolemtyCoordinates) for example.

=== Arguments ===

M — manifold to check this for epsilon — maximal allowed error when checking the relations, use None for exact comparison.

evaluate_word(word, G=None)

Given a word in the generators of the fundamental group, compute the corresponding matrix. By default, these are the generators of the unsimplified presentation of the fundamental group. An optional SnapPy fundamental group can be given if the words are in generators of a different presentation, e.g., c.evaluate_word(word, M.fundamental_group(True)) to evaluate a word in the simplified presentation returned by M.fundamental_group(True).

For now, the matrix is returned as list of lists.

static from_snappy_manifold(M, dec_prec=None, bits_prec=None, intervals=False)

Constructs an assignment of shape parameters/cross ratios using the tetrahehdra_shapes method of a given SnapPy manifold. The optional parameters are the same as that of tetrahedra_shapes.

get_manifold()

Get the manifold for which this structure represents a solution to the gluing equations.

induced_representation(N)

Given a PSL(2,C) representation constructs the induced representation for the given N. The induced representation is in SL(N,C) if N is odd and SL(N,C) / {+1,-1} if N is even and is described in the Introduction of Garoufalidis, Thurston, Zickert The Complex Volume of SL(n,C)-Representations of 3-Manifolds http://arxiv.org/abs/1111.2828

There is a canonical group homomorphism SL(2,C)->SL(N,C) coming from the the natural SL(2,C)-action on the vector space Sym^{N-1}(C^2). This homomorphisms decends to a homomorphism from PSL(2,C) if one divides the right side by {+1,-1} when N is even. Composing a representation with this homomorphism gives the induced representation.

is_geometric(epsilon=1e-06)

Returns true if all shapes corresponding to this solution have positive imaginary part.

If the solutions are exact, it returns true if one of the corresponding numerical solutions is geometric.

An optional epsilon can be given. An imaginary part of a shape is considered positive if it is larger than this epsilon.

is_induced_from_psl2(epsilon=None)

For each simplex and each edges, checks that all cross ratios of that simplex that are parallel to that each are the same (maximal absolute difference is the epsilon given as argument). This means that the corresponding representation is induced by a PSL(2,C) representation.

is_pu_2_1_representation(epsilon, epsilon2=None)

Returns True if the representation is also a PU(2,1)-representation. This uses Proposition 3.5 and the remark following that proposition in [FKR2013].

If a condition given in that Proposition is violated, the method returns an object whose Boolean value is still False and that indicates which condition was violated. Thus, this method can still be used in if statements.

The method tests the following complex equalities and inequalities:

  • the three complex equations given in (3.5.1) of [FKR2013].

  • the inequality zijl \(\\not=\) -1.

Remark: It does not check whether all zij * zji are real or not as these are still valid CR configurations, see the remark following Proposition 3.5.

The user has to supply an epsilon: an equality/inequality is considered to be true if and only if the absolute value | LHS - RHS | of difference between the left and right hand side is less/greater than epsilon.

The user can supply another parameter, epsilon2. If any | LHS - RHS | is in the interval [epsilon, epsilon2], this method fails with an exception as the value of | LHS - RHS | is an ambiguous interval where it is unclear whether inequality fails to hold because it truly does hold or just because of numerical noise.

is_real(epsilon)

Returns True if all cross ratios are real (have absolute imaginary part < epsilon where epsilon is given as argument). This means that the corresponding representation is in PSL(N,R).

long_edge(tet, v0, v1, v2)

The matrix that labels a long edge starting at vertex (v0, v1, v2) of a doubly truncated simplex corresponding to the ideal tetrahedron with index tet.

This matrix was labeled alpha^{v0v1v2} in Figure 18 of Garoufalidis, Goerner, Zickert: Gluing Equations for PGL(n,C)-Representations of 3-Manifolds http://arxiv.org/abs/1207.6711

It is computed using equation 10.22.

The resulting matrix is given as a python list of lists.

middle_edge(tet, v0, v1, v2)

The matrix that labels a middle edge starting at vertex (v0, v1, v2) of a doubly truncated simplex corresponding to the ideal tetrahedron with index tet.

This matrix was labeled beta^{v0v1v2} in Figure 18 of Garoufalidis, Goerner, Zickert: Gluing Equations for PGL(n,C)-Representations of 3-Manifolds http://arxiv.org/abs/1207.6711

It is computed using equation 10.22.

The resulting matrix is given as a python list of lists.

multiply_and_simplify_terms_in_RUR()

If a cross ratio is given as Rational Univariate Representation with numerator and denominator being a product, multiply the terms, reduce the fraction and return the result.

See multiply_and_simplify_terms of RUR.

This loses information about how the numerator and denominator are factorised.

multiply_terms_in_RUR()

If a cross ratio is given as Rational Univariate Representation with numerator and denominator being a product, multiply the terms and return the result.

See multiply_terms of RUR.

This loses information about how the numerator and denominator are factorised.

num_tetrahedra()

The number of tetrahedra for which we have cross ratios.

numerical()

Turn exact solutions into numerical solutions using pari. Similar to numerical() of PtolemyCoordinates. See help(ptolemy.PtolemyCoordinates) for example.

short_edge(tet, v0, v1, v2)

The matrix that labels a long edge starting at vertex (v0, v1, v2) of a doubly truncated simplex corresponding to the ideal tetrahedron with index tet.

This matrix was labeled gamma^{v0v1v2} in Figure 18 of Garoufalidis, Goerner, Zickert: Gluing Equations for PGL(n,C)-Representations of 3-Manifolds http://arxiv.org/abs/1207.6711

It is computed using equation 10.22.

The resulting matrix is given as a python list of lists.

to_PUR()

If any Ptolemy coordinates are given as Rational Univariate Representation, convert them to Polynomial Univariate Representation and return the result.

See to_PUR of RUR.

This conversion might lead to very large coefficients.

volume_numerical(drop_negative_vols=False)

Turn into (Galois conjugate) numerical solutions and compute volumes. If already numerical, only compute the one volume. See numerical().

If drop_negative_vols = True is given as optional argument, only return non-negative volumes.

x_coordinate(tet, pt)

Returns the X-coordinate for the tetrahedron with index tet at the point pt (quadruple of integers adding up to N).

See Definition 10.9: Garoufalidis, Goerner, Zickert: Gluing Equations for PGL(n,C)-Representations of 3-Manifolds http://arxiv.org/abs/1207.6711

Flattenings

class snappy.ptolemy.coordinates.Flattenings(d, manifold_thunk=<function Flattenings.<lambda>>, evenN=2)

Represents a flattening assigned to each edge of a simplex as dictionary.

We assign to each pair of parallel edges of each simplex a triple (w, z, p) such that:

w = log(z) + p * (2 * pi * i / N)   where N is fixed and even.

For N = 2, the three triples belonging to a simplex form a combinatorial flattening (w0, w1, w2) as defined in Definition 3.1 in Walter D. Neumann, Extended Bloch group and the Cheeger-Chern-Simons class http://arxiv.org/abs/math.GT/0307092

For N > 2, the three triples form a generalized combinatorial flattening (w0, w1, w2) that gives an element in the generalized Extended Bloch group which is the Extended Bloch group corresponding to the Riemann surface given by:

u1 * e^w0 + u2 * e^w1 = 1

where u1^N = u2^N = 1.

A representation in SL(n,C) and SL(n,C)/{+1,-1} with n even gives an element in the generalized Extended Bloch group for N = 2. A representation in PSL(n,C) with n even in the group for N = n. A representation in PSL(n,C) with n odd in the group for N = 2 * n.

This work has not been published yet.

If f is a flattening, then in the notation of Neumann, the value of:

f['z_xxxx_y']    is (w0, z, p)
f['zp_xxxx_y']   is (w1, z', q)
f['zpp_xxxx_y']  is (w2, z'', r).
N()

Get the N such that these cross ratios are for SL/PSL(N,C)-representations.

check_against_manifold(M=None, epsilon=1e-10)

Checks that the flattening really is a solution to the logarithmic PGL(N,C) gluing equations of a manifold. Usage similar to check_against_manifold of Ptolemy Coordinates, see help(ptolemy.Coordinates) for similar examples.

=== Arguments ===

M — manifold to check this for epsilon — maximal allowed error when checking the equations

complex_volume(with_modulo=False)

Compute complex volume. The complex volume is defined only up to some multiple of m where m = i * pi**2/6 for PSL(2,C) and SL(N,C) and m = i * pi**2/18 for PSL(3,C).

When called with with_modulo = True, gives a pair (volume, m)

classmethod from_tetrahedra_shapes_of_manifold(M)

Takes as argument a manifold and produces (weak) flattenings using the tetrahedra_shapes of the manifold M.

>>> from snappy import Manifold
>>> M = Manifold("5_2")
>>> flattenings = Flattenings.from_tetrahedra_shapes_of_manifold(M)
>>> flattenings.check_against_manifold(M)
>>> flattenings.check_against_manifold()
get_manifold()

Get the manifold for which this structure represents a flattening.

get_order()

Returns the number N. This flattening represents an element in the generalized Extended Bloch group for the Riemann surface given by u1 * e^w0 + u2 * e^w1 = 1 where u1^N = u2^N = 1.

get_zpq_triple(key_z)

Gives a flattening as triple [z;p,q] representing an element in the generalized Extended Bloch group similar to the way the triple [z;p,q] is used in Lemma 3.2 in Walter D. Neumann, Extended Bloch group and the Cheeger-Chern-Simons class http://arxiv.org/abs/math.GT/0307092

num_tetrahedra()

The number of tetrahedra for which we have cross ratios.

NonZeroDimensionalComponent

class snappy.ptolemy.component.NonZeroDimensionalComponent(witnesses=[], dimension='unknown', free_variables=None, genus=None, p=None)

Represents a non-zero dimensional component in the Ptolemy variety. It is a list that can hold points sampled from that component (witnesses).

Other functions

snappy.ptolemy.solutions_from_magma(output, numerical=False)

Obsolete, use processFileDispatch.parse_solutions instead.

Assumes the given string is the output of a magma computation, parses it and returns a list of solutions. A non-zero dimensional component of the variety is reported as NonZeroDimensionalComponent.

snappy.ptolemy.solutions_from_magma_file(filename, numerical=False)

Obsolete, use processFileDispatch.parse_solutions_from_file instead.

Reads the output from a magma computation from the file with the given filename and returns a list of solutions. Also see solutions_from_magma. A non-zero dimensional component of the variety is reported as NonZeroDimensionalComponent.

[FKR2013] (1,2)

Falbel, Koseleff, Rouillier: Representations of Fundamental Groups of 3-Manifolds into PGL(3,C): Exact Computations in Low Complexity, http://arxiv.org/abs/1307.6697