Internals of verified computations
Naming
The names of methods containing check
will raise an exception if
the desired property cannot be certified. There are different types of
Exceptions to indicate how the certification failed. This type can be
used by other methods to perform some action such as changing the
triangulation or increasing precision or to give up.
The user-facing methods have names starting with verify
or
verified
and will fail more gracefully returning False
or
None
in such a case.
Generating certified shape intervals
The recommended way to obtain certified intervals for the shapes is via
manifold.tetrahedra_shapes(intervals=True)
as described
earlier. Here we document the KrawczykShapesEngine
and
IntervalNewtonShapesEngine
which is implemented internally to
generate the intervals. It is of interest for those users who want to
understand the underlying interval math and experiment with the Newton
interval method or the Krawczyk test. CertifiedShapesEngine
is an
alias of either KrawczykShapesEngine
or
IntervalNewtonShapesEngine
to determine the default method used by
verify.
- snappy.verify.CertifiedShapesEngine
alias of
KrawczykShapesEngine
- class snappy.verify.IntervalNewtonShapesEngine(M, initial_shapes, bits_prec=None, dec_prec=None)
An engine that is initialized with an approximated candidate solution to the rectangular gluing equations and produces intervals certified to contain a true solution. After the engine is successfully run, the resulting intervals are stored in certified_shapes which is a vector of elements in a Sage’s ComplexIntervalField.
A simple example to obtain certified shape intervals that uses KrawczykShapesEngine or IntervalNewtonShapesEngine under the hood:
sage: from snappy import Manifold sage: M = Manifold("m015") sage: M.tetrahedra_shapes('rect', bits_prec = 80, intervals = True) # doctest: +NUMERIC15 +NORMALIZE_WHITESPACE [0.6623589786223730129805? + 0.5622795120623012438992?*I, 0.6623589786223730129805? + 0.5622795120623012438992?*I, 0.6623589786223730129805? + 0.5622795120623012438992?*I]
Its objective is thus the same as HIKMOT and it is certainly HIKMOT inspired. However, it conceptually differs in that:
It uses the Newton interval method instead of the Krawczyk test (we implement Gaussian elimination in interval arithmetic to compute the inverse of an interval matrix having interval arithmetic semantics, see mat_solve).
It uses complex numbers in it’s Newton interval method. We simply use Sage’s complex interval type avoiding the need of converting n x n complex matrices into 2n x 2n real matrices as described Section 3.4 of the HIKMOT paper.
We avoid automatic differentiation. We pick an independent set of equations of the following form and try to solve them:
log(LHS) = 0
where
LHS = c * z0^a0 * (1-z0)^b0 * z1^a1 * (1-z1)^b1 * …
with a, b and c’s as returned by Manifold.gluing_equations(‘rect’).
The derivative of log (LHS) with respect to zj is simply given by
aj/zj - bj/(1-zj)
and thus no need for automatic differentiation.
In contrast to HIKMOT, we use and return Sage’s native implementation of (complex) interval arithmetic here, which allows for increased interoperability. Another advantage is that Sage supports arbitrary precision. Unfortunately, performance suffers and this implementation is 5-10 times slower than HIKMOT.
Here is an example how to explicitly invoke the IntervalNewtonShapesEngine:
sage: shapes = M.tetrahedra_shapes('rect', bits_prec = 80) sage: C = IntervalNewtonShapesEngine(M, shapes, bits_prec = 80) sage: C.expand_until_certified() True sage: C.certified_shapes # doctest: +ELLIPSIS (0.662358978622373012981? + 0.562279512062301243...?*I, 0.66235897862237301298...? + 0.562279512062301243...?*I, 0.66235897862237301298...? + 0.562279512062301243...?*I)
- static certified_newton_iteration(equations, shape_intervals, point_in_intervals=None, interval_value_at_point=None)
Given shape intervals z, performs a Newton interval iteration N(z) as described in newton_iteration. Returns a pair (boolean, N(z)) where the boolean is True if N(z) is contained in z.
If the boolean is True, it is certified that N(z) contains a true solution, e.g., a point for which f is truly zero.
See newton_iteration for the other parameters.
This follows from Theorem 1 of Zgliczynski’s notes.
Some examples:
sage: from snappy import Manifold sage: M = Manifold("m019") sage: C = IntervalNewtonShapesEngine(M, M.tetrahedra_shapes('rect'), ... bits_prec = 80)
Intervals containing the true solution:
sage: from sage.all import vector sage: good_shapes = vector([ ... C.CIF(C.RIF(0.78055, 0.78056), C.RIF(0.91447, 0.91448)), ... C.CIF(C.RIF(0.78055, 0.78056), C.RIF(0.91447, 0.91448)), ... C.CIF(C.RIF(0.46002, 0.46003), C.RIF(0.63262, 0.63263))]) sage: is_certified, shapes = IntervalNewtonShapesEngine.certified_newton_iteration(C.equations, good_shapes) sage: is_certified True sage: shapes # doctest: +ELLIPSIS (0.78055253? + 0.91447366...?*I, 0.7805525...? + 0.9144736...?*I, 0.4600211...? + 0.632624...?*I)
This means that a true solution to the rectangular gluing equations is contained in both the given intervals (good_shapes) and the returned intervals (shapes) which are a refinement of the given intervals.
Intervals not containing a true solution:
sage: from sage.all import vector sage: bad_shapes = vector([ ... C.CIF(C.RIF(0.78054, 0.78055), C.RIF(0.91447, 0.91448)), ... C.CIF(C.RIF(0.78055, 0.78056), C.RIF(0.91447, 0.91448)), ... C.CIF(C.RIF(0.46002, 0.46003), C.RIF(0.63262, 0.63263))]) sage: is_certified, shapes = IntervalNewtonShapesEngine.certified_newton_iteration(C.equations, bad_shapes) sage: is_certified False
- expand_until_certified(verbose=False)
Try Newton interval iterations, expanding the shape intervals until we can certify they contain a true solution. If succeeded, return True and write certified shapes to certified_shapes. Set verbose = True for printing additional information.
- static interval_vector_is_contained_in(vecA, vecB)
Given two vectors of intervals, return whether the first one is contained in the second one. Examples:
sage: RIF = RealIntervalField(80) sage: CIF = ComplexIntervalField(80) sage: box = CIF(RIF(-1,1),RIF(-1,1)) sage: a = [ CIF(0.1), CIF(1) + box ] sage: b = [ CIF(0) + box, CIF(1) + 2 * box ] sage: c = [ CIF(0), CIF(1) + 3 * box ] sage: IntervalNewtonShapesEngine.interval_vector_is_contained_in(a, b) True sage: IntervalNewtonShapesEngine.interval_vector_is_contained_in(a, c) False sage: IntervalNewtonShapesEngine.interval_vector_is_contained_in(b, a) False sage: IntervalNewtonShapesEngine.interval_vector_is_contained_in(b, c) False sage: IntervalNewtonShapesEngine.interval_vector_is_contained_in(c, a) False sage: IntervalNewtonShapesEngine.interval_vector_is_contained_in(c, b) False
- static interval_vector_mid_points(vec)
Given a vector of complex intervals, return the midpoints (as 0-length complex intervals) of them.
- static interval_vector_union(vecA, vecB)
Given two vectors of intervals, return the vector of their unions, i.e., the smallest interval containing both intervals.
- static log_gluing_LHS_derivatives(equations, shapes)
Compute the Jacobian of the vector-valued function f described in the above log_gluing_LHSs:
sage: from snappy import Manifold sage: M = Manifold("m019") sage: equations = M.gluing_equations('rect') sage: RIF = RealIntervalField(80) sage: CIF = ComplexIntervalField(80) sage: shape1 = CIF(RIF(0.78055,0.78056), RIF(0.9144, 0.9145)) sage: shape2 = CIF(RIF(0.46002,0.46003), RIF(0.6326, 0.6327)) sage: shapes = [shape1, shape1, shape2] sage: IntervalNewtonShapesEngine.log_gluing_LHS_derivatives(equations, shapes) # doctest: +ELLIPSIS +NORMALIZE_WHITESPACE [ 0.292? - 1.66...?*I 0.292? - 1.66...?*I 0.752? - 1.034...?*I] [-0.5400? + 0.63...?*I -0.5400? + 0.63...?*I 1.561? + 1.829...?*I] [ 0.2482? + 1.034...?*I 0.2482? + 1.034...?*I -2.313? - 0.795...?*I] [ 0.5400? - 0.63...?*I -0.5400? + 0.63...?*I 0] [...-0.4963? - 2.068?*I 1.0800? - 1.26...?*I 0.752? - 1.034...?*I]
- static log_gluing_LHSs(equations, shapes)
Given the result of M.gluing_equations(‘rect’) or a subset of rows of it and shapes, return a vector of log(LHS) where
LHS = c * z0 ** a0 * (1-z0) ** b0 * z1 ** a1 * …
Let f: C^n -> C^n denote the function which takes shapes and returns the vector of log(LHS).
The reason we take the logarithm of the rectangular gluing equations is because the logarithmic derivative is of a particular nice form:
sage: from snappy import Manifold sage: M = Manifold("m019") sage: equations = M.gluing_equations('rect') sage: RIF = RealIntervalField(80) sage: CIF = ComplexIntervalField(80) sage: zero = CIF(0).center() sage: shape1 = CIF(RIF(0.78055,0.78056), RIF(0.9144, 0.9145)) sage: shape2 = CIF(RIF(0.46002,0.46003), RIF(0.6326, 0.6327))
An interval solution containing the true solution. The log of each rectangular equation should be 0 for the true solution, hence the interval should contain zero:
sage: shapes = [shape1, shape1, shape2] sage: LHSs = IntervalNewtonShapesEngine.log_gluing_LHSs(equations, shapes) sage: LHSs # doctest: +ELLIPSIS (0.000? + 0.000?*I, 0.000? + 0.000?*I, 0.000? + 0.000?*I, 0.000...? + 0.000...?*I, 0.000? + 0.000?*I) sage: zero in LHSs[0] True
An interval not containing the true solution:
sage: shapes = [shape1, shape1, shape1] sage: LHSs = IntervalNewtonShapesEngine.log_gluing_LHSs(equations, shapes) sage: LHSs # doctest: +ELLIPSIS (0.430? - 0.078?*I, -0.2...? + 0.942?*I, -0.1...? - 0.8...?*I, 0.000...? + 0.000...?*I, 0.430? - 0.078?*I) sage: zero in LHSs[0] False
- static newton_iteration(equations, shape_intervals, point_in_intervals=None, interval_value_at_point=None)
Perform a Newton interval method of iteration for the function f described in log_gluing_LHSs.
Let z denote the shape intervals. Let z_center be a point close to the center point of the shape intervals (in the implementation, z_center is an interval of again, of length zero).
The result returned will be
N(z) = z_center - ((Df)(z))^-1 f(z_center)
The user can overwrite the z_center to be used by providing point_in_intervals (which have to be 0-length complex intervals). The user can also give the interval value of f(z_center) by providing interval_value_at_point to avoid re-evaluation of f(z_center).
A very approximate solution:
sage: from snappy import Manifold sage: M = Manifold("m019") sage: shapes = [ 0.7+1j, 0.7+1j, 0.5+0.5j ]
Get the equations and initialize zero-length intervals from it:
sage: C = IntervalNewtonShapesEngine(M, shapes, bits_prec = 80) sage: C.initial_shapes (0.69999999999999995559107902? + 1*I, 0.69999999999999995559107902? + 1*I, 0.50000000000000000000000000? + 0.50000000000000000000000000?*I)
Do several Newton interval operations to get a better solution:
sage: shape_intervals = C.initial_shapes sage: for i in range(4): # doctest: +ELLIPSIS ... shape_intervals = IntervalNewtonShapesEngine.newton_iteration(C.equations, shape_intervals) ... print(shape_intervals) (0.78674683118381457770...? + 0.9208680745160821379529?*I, 0.786746831183814577703...? + 0.9208680745160821379529?*I, 0.459868058287098030934...? + 0.61940871855835167317...?*I) (0.78056102517632648594...? + 0.9144962118446750482...?*I, 0.78056102517632648594...? + 0.9144962118446750482...?*I, 0.4599773577869384936554? + 0.63251940718694538695...?*I) (0.78055253104531610049...? + 0.9144736621585220345231?*I, 0.780552531045316100497...? + 0.9144736621585220345231?*I, 0.460021167103732494700...? + 0.6326241909236695020810...?*I) (0.78055252785072483256...? + 0.91447366296772644033...?*I, 0.7805525278507248325678? + 0.914473662967726440333...?*I, 0.4600211755737178641204...? + 0.6326241936052562241142...?*I)
For comparison:
sage: M.tetrahedra_shapes('rect') [0.780552527850725 + 0.914473662967726*I, 0.780552527850725 + 0.914473662967726*I, 0.460021175573718 + 0.632624193605256*I]
Start with a rather big interval, note that the Newton interval method is stable in the sense that the interval size decreases:
sage: box = C.CIF(C.RIF(-0.0001,0.0001),C.RIF(-0.0001,0.0001)) sage: shape_intervals = C.initial_shapes.apply_map(lambda shape: shape + box) sage: shape_intervals (0.700? + 1.000?*I, 0.700? + 1.000?*I, 0.500? + 0.500?*I) sage: for i in range(7): ... shape_intervals = IntervalNewtonShapesEngine.newton_iteration(C.equations, shape_intervals) sage: print(shape_intervals) # doctest: +ELLIPSIS (0.78055252785072483798...? + 0.91447366296772645593...?*I, 0.7805525278507248379869? + 0.914473662967726455938...?*I, 0.460021175573717872891...? + 0.632624193605256171637...?*I)
- class snappy.verify.KrawczykShapesEngine(M, initial_shapes, bits_prec=None, dec_prec=None)
An engine that is initialized with an approximated candidate solution to the rectangular gluing equations and produces intervals certified to contain a true solution. After the engine is successfully run, the resulting intervals are stored in certified_shapes which is a vector of elements in a Sage’s ComplexIntervalField.
A simple example to obtain certified shape intervals that uses the KrawczykShapesEngine or IntervalNewtonShapesEngine under the hood:
sage: from snappy import Manifold sage: M = Manifold("m015") sage: M.tetrahedra_shapes('rect', bits_prec = 80, intervals = True) # doctest: +NUMERIC15 +NORMALIZE_WHITESPACE [0.6623589786223730129805? + 0.5622795120623012438992?*I, 0.6623589786223730129805? + 0.5622795120623012438992?*I, 0.6623589786223730129805? + 0.5622795120623012438992?*I]
Its objective is thus the same as HIKMOT and it is certainly HIKMOT inspired. However, it conceptually differs in that:
It uses complex numbers in it’s computations. We simply use Sage’s complex interval type avoiding the need of converting n x n complex matrices into 2n x 2n real matrices as described Section 3.4 of the HIKMOT paper.
We avoid automatic differentiation. We pick an independent set of equations of the following form and try to solve them:
log(LHS) = 0
where
LHS = c * z0^a0 * (1-z0)^b0 * z1^a1 * (1-z1)^b1 * …
with a, b and c’s as returned by Manifold.gluing_equations(‘rect’).
The derivative of log (LHS) with respect to zj is simply given by
aj/zj - bj/(1-zj)
and thus no need for automatic differentiation.
For speed-up, the approximate inverse is always computed with double’s. Some intermediate matrix computations are performed sparsely.
In contrast to HIKMOT, we use and return Sage’s native implementation of (complex) interval arithmetic here, which allows for increased interoperability. Another advantage is that Sage supports arbitrary precision.
Here is an example how to explicitly invoke the KrawczykShapesEngine:
sage: shapes = M.tetrahedra_shapes('rect', bits_prec = 80) sage: C = KrawczykShapesEngine(M, shapes, bits_prec = 80) sage: C.expand_until_certified() True sage: C.certified_shapes # doctest: +NUMERIC12 (0.6623589786223730129805? + 0.5622795120623012438992?*I, 0.6623589786223730129805? + 0.5622795120623012438992?*I, 0.6623589786223730129805? + 0.5622795120623012438992?*I)
And here an example where the initial solution is somewhat off:
sage: M = Manifold("m019") sage: shapes = [0.780559+0.91449j, 0.780559+0.9144j, 0.46009+0.639j] sage: C = KrawczykShapesEngine(M, shapes, bits_prec = 100) sage: C.expand_until_certified() True sage: C.certified_shapes # doctest: +ELLIPSIS (0.7806? + 0.9145?*I, 0.7806? + 0.9145?*I, 0.460...? + 0.6326?*I)
- expand_until_certified(verbose=False)
Try Krawczyk iterations (i.e., expanding the shape intervals [z] by the Krawczyk interval K(z0, [z], f)) until we can certify they contain a true solution.
If succeeded, return True and write certified shapes to certified_shapes. Set verbose = True for printing additional information.
- static interval_vector_is_contained_in(vecA, vecB)
Given two vectors of intervals, return whether the first one is contained in the second one. Examples:
sage: RIF = RealIntervalField(80) sage: CIF = ComplexIntervalField(80) sage: box = CIF(RIF(-1,1),RIF(-1,1)) sage: a = [ CIF(0.1), CIF(1) + box ] sage: b = [ CIF(0) + box, CIF(1) + 2 * box ] sage: c = [ CIF(0), CIF(1) + 3 * box ] sage: KrawczykShapesEngine.interval_vector_is_contained_in(a, b) True sage: KrawczykShapesEngine.interval_vector_is_contained_in(a, c) False sage: KrawczykShapesEngine.interval_vector_is_contained_in(b, a) False sage: KrawczykShapesEngine.interval_vector_is_contained_in(b, c) False sage: KrawczykShapesEngine.interval_vector_is_contained_in(c, a) False sage: KrawczykShapesEngine.interval_vector_is_contained_in(c, b) False
- static interval_vector_mid_points(vec)
Given a vector of complex intervals, return the midpoints (as 0-length complex intervals) of them.
- static interval_vector_union(vecA, vecB)
Given two vectors of intervals, return the vector of their unions, i.e., the smallest interval containing both intervals.
- krawczyk_interval(shape_intervals)
Compute the interval in the Krawczyk test.
It is given as
K(z0, [z], f) := z0 - c * f(z0) + (Id - c * df([z])) * ([z] - z0)
- where
z0 is the approximate candidate solution,
[z] are the shape_intervals we try to verify,
f is the function taking the shapes to the errors of the logarithmic gluing equations
c is an approximate inverse of df
df([z]) is the derivative of f (interval-)evaluated for [z]
Note that z0 in self.initial_shapes which are complex intervals containing only one value (the candidate solution given initially).
If K is contained in [z], then we have proven that [z] contains a solution to the gluing equations.
Do several Krawczyk operations to get a better solution:
sage: from sage.all import vector sage: M = Manifold("m019") sage: shapes = vector(ComplexIntervalField(53), [ 0.5+0.8j, 0.5+0.8j, 0.5+0.8j]) sage: for i in range(15): ... penultimateShapes = shapes ... centers = [ shape.center() for shape in shapes ] ... C = KrawczykShapesEngine(M, centers, bits_prec = 53) ... shapes = C.krawczyk_interval(shapes) sage: shapes # doctest: +NUMERIC12 (0.78055252785073? + 0.91447366296773?*I, 0.780552527850725? + 0.91447366296773?*I, 0.460021175573718? + 0.632624193605256?*I)
- log_gluing_LHS_derivatives(shapes)
Compute the Jacobian of the vector-valued function f described in the above log_gluing_LHSs:
sage: from snappy import Manifold sage: M = Manifold("m019") sage: shapes = M.tetrahedra_shapes('rect', bits_prec = 80) sage: C = KrawczykShapesEngine(M, shapes, bits_prec = 80) sage: RIF = RealIntervalField(80) sage: CIF = ComplexIntervalField(80) sage: shape1 = CIF(RIF(0.78055,0.78056), RIF(0.9144, 0.9145)) sage: shape2 = CIF(RIF(0.46002,0.46003), RIF(0.6326, 0.6327)) sage: shapes = [shape1, shape1, shape2] sage: C.log_gluing_LHS_derivatives(shapes) # doctest: +NUMERIC3 [ 0.292? - 1.6666?*I 0.292? - 1.6666?*I 0.752? - 1.0340?*I] [ 0.5400? - 0.6327?*I 0.5400? - 0.6327?*I -1.561? - 1.8290?*I] [ 0.5400? - 0.6327?*I -0.5400? + 0.6327?*I 0]
- log_gluing_LHS_derivatives_sparse(shapes)
A column-sparse matrix version of log_gluing_LHS_derivatives_sparse. The result is a list of list of pairs. Each list of pairs corresponds to a column, a pair being (index of row, value) where the index is increasing.
- log_gluing_LHSs(shapes)
Given the result of M.gluing_equations(‘rect’) or a subset of rows of it and shapes, return a vector of log(LHS) where
LHS = c * z0 ** a0 * (1-z0) ** b0 * z1 ** a1 * …
Let f: C^n -> C^n denote the function which takes shapes and returns the vector of log(LHS).
The reason we take the logarithm of the rectangular gluing equations is because the logarithmic derivative is of a particular nice form:
sage: from snappy import Manifold sage: M = Manifold("m019") sage: equations = M.gluing_equations('rect') sage: RIF = RealIntervalField(80) sage: CIF = ComplexIntervalField(80) sage: zero = CIF(0).center() sage: shape1 = CIF(RIF(0.78055,0.78056), RIF(0.9144, 0.9145)) sage: shape2 = CIF(RIF(0.46002,0.46003), RIF(0.6326, 0.6327))
An interval solution containing the true solution. The log of each rectangular equation should be 0 for the true solution, hence the interval should contain zero:
sage: shapes = [shape1, shape1, shape2] sage: C = KrawczykShapesEngine(M, [shape.center() for shape in shapes], bits_prec = 53) sage: LHSs = C.log_gluing_LHSs(shapes) sage: LHSs # doctest: +NUMERIC6 (0.000? + 0.000?*I, 0.000? + 0.000?*I, 0.0000? + 0.0000?*I) sage: zero in LHSs[0] True
An interval not containing the true solution:
sage: shapes = [shape1, shape1, shape1] sage: LHSs = C.log_gluing_LHSs(shapes) sage: LHSs # doctest: +NUMERIC3 (0.430? - 0.078?*I, 0.246? - 0.942?*I, 0.0000? + 0.0000?*I) sage: zero in LHSs[0] False
- static matrix_times_sparse(m, sparse_m)
Multiply a (dense) Sage matrix with a column-sparse matrix (in the format described in log_gluing_LHS_derivatives_sparse).
Verification of hyperbolicity
Methods containing check
will raise an exception if the desired property
cannot be certified. Methods containing verify
or verified
will fail
more gracefully returning False
or None
in such a case.
- snappy.verify.hyperbolicity.check_logarithmic_gluing_equations_and_positively_oriented_tets(manifold, shape_intervals)
Given a SnapPy manifold manifold and complex intervals for the shapes shape_intervals that are certified to contain a solution to the rectangular gluing equations, verify that the logarithmic gluing equations are also fulfilled and that all shapes have positive imaginary part. It will raise an exception if the verification fails. This is sufficient to prove that the manifold is indeed hyperbolic.
Since the given interval are supposed to contain a true solution of the rectangular gluing equations, the logarithmic gluing equations are known to be fulfilled up to a multiple of 2 pi i. Thus it is enough to certify that the absolute error of the logarithmic gluing equations is < 0.1. Using interval arithmetic, this function certifies this and positivity of the imaginary parts of the shapes:
sage: from snappy import Manifold sage: M = Manifold("m019") sage: check_logarithmic_gluing_equations_and_positively_oriented_tets( ... M, M.tetrahedra_shapes('rect', intervals=True))
The SnapPy triangulation of the following hyperbolic manifold contains actually negatively oriented tetrahedra:
sage: M = Manifold("t02774") sage: check_logarithmic_gluing_equations_and_positively_oriented_tets( ... M, M.tetrahedra_shapes('rect', intervals=True)) # doctest: +IGNORE_EXCEPTION_DETAIL Traceback (most recent call last): ... ShapePositiveImaginaryPartNumericalVerifyError: Numerical verification that shape has positive imaginary part has failed: Im(0.4800996900657? - 0.0019533695046?*I) > 0
Verified canonical cell decompositions
- snappy.verify.canonical.interval_checked_canonical_triangulation(M, bits_prec=None)
Given a canonical triangulation of a cusped (possibly non-orientable) manifold M, return this triangulation if it has tetrahedral cells and can be verified using interval arithmetics with the optional, given precision. Otherwise, raises an Exception.
It fails when we call it on something which is not the canonical triangulation:
sage: from snappy import Manifold sage: M = Manifold("m015") sage: interval_checked_canonical_triangulation(M) # doctest: +ELLIPSIS +IGNORE_EXCEPTION_DETAIL Traceback (most recent call last): ... TiltProvenPositiveNumericalVerifyError: Numerical verification that tilt is negative has failed, tilt is actually positive. This is provably not the proto-canonical triangulation: 0.164542163...? <= 0
It verifies the canonical triangulation:
sage: M.canonize() sage: K = interval_checked_canonical_triangulation(M) sage: K m015(0,0)
Has a non-tetrahedral canonical cell:
sage: M = Manifold("m137") sage: M.canonize() sage: interval_checked_canonical_triangulation(M) # doctest: +ELLIPSIS +IGNORE_EXCEPTION_DETAIL Traceback (most recent call last): ... TiltInequalityNumericalVerifyError: Numerical verification that tilt is negative has failed: 0.?e-1... < 0
Has a cubical canonical cell:
sage: M = Manifold("m412") sage: M.canonize() sage: interval_checked_canonical_triangulation(M) # doctest: +ELLIPSIS +IGNORE_EXCEPTION_DETAIL Traceback (most recent call last): ... TiltInequalityNumericalVerifyError: Numerical verification that tilt is negative has failed: 0.?e-1... < 0
- snappy.verify.canonical.exactly_checked_canonical_retriangulation(M, bits_prec, degree)
Given a proto-canonical triangulation of a cusped (possibly non-orientable) manifold M, return its canonical retriangulation which is computed from exact shapes. The exact shapes are computed using snap (which uses the LLL-algorithm). The precision (in bits) and the maximal degree need to be specified (here 300 bits precision and polynomials of degree less than 4):
sage: from snappy import Manifold sage: M = Manifold("m412") sage: M.canonize() sage: K = exactly_checked_canonical_retriangulation(M, 300, 4)
M’s canonical cell decomposition has a cube, so non-tetrahedral:
sage: K.has_finite_vertices() True
Has 12 tetrahedra after the retrianglation:
sage: K.num_tetrahedra() 12
Check that it fails on something which is not a proto-canonical triangulation:
sage: from snappy import Manifold sage: M = Manifold("m015") sage: exactly_checked_canonical_retriangulation(M, 500, 6) # doctest: +IGNORE_EXCEPTION_DETAIL Traceback (most recent call last): ... TiltProvenPositiveNumericalVerifyError: Numerical verification that tilt is negative has failed, tilt is actually positive. This is provably not the proto-canonical triangulation: 0.1645421638874662848910671879? <= 0
Exact computations for cusp cross sections
The square_extensions module provides two special classes to give exact representations of the values involved when computing a cusp cross section.
The method find_shapes_as_complex_sqrt_lin_combinations returns a list of shapes as ComplexSqrtLinCombination’s. This can be used as input to CuspCrossSection. The outputs of CuspCrossSection, including the tilts, will then be of type SqrtLinCombination.
Consider the real number field N generated by the real and imaginary part of the shapes. The edge lengths and the factors used to normalize the cusp areas will be square roots in N and thus the tilts will be N-linear combinations of square roots in N. To avoid computing in a massive tower of square extensions of N, we implement SqrtLinCombination here that provides a special implementation of the == operator.
- snappy.verify.square_extensions.find_shapes_as_complex_sqrt_lin_combinations(M, prec, degree)
Given a manifold M, use snap (which uses LLL-algorithm) with the given decimal precision and maximal degree to find exact values for the shapes’ real and imaginary part. Return the shapes as list of ComplexSqrtLinCombination’s. Return None on failure.
Example:
sage: from snappy import Manifold sage: M=Manifold("m412") sage: find_shapes_as_complex_sqrt_lin_combinations(M, 200, 10) [ComplexSqrtLinCombination((1/2) * sqrt(1), (x - 1/2) * sqrt(1)), ComplexSqrtLinCombination((1/2) * sqrt(1), (x - 1/2) * sqrt(1)), ComplexSqrtLinCombination((1/2) * sqrt(1), (x - 1/2) * sqrt(1)), ComplexSqrtLinCombination((1/2) * sqrt(1), (x - 1/2) * sqrt(1)), ComplexSqrtLinCombination((1/2) * sqrt(1), (x - 1/2) * sqrt(1))]
- class snappy.verify.square_extensions.SqrtLinCombination(value=None, d={}, embed_cache=None)
A class representing a linear combination
c_1 * sqrt(r_1) + c_2 * sqrt(r_2) + … + c_n * sqrt(r_n)
where c_i and r_i have to be of type Integer, Rational or elements of the same Sage NumberField with a real embedding (Caution: this is assumed but not checked!) such that all r_i are positive (Caution: this is not checked during construction!).
It implements +, -, * where one of the operators is allowed to be an integer or rational.
/ is only implemented when the denominator has only one term c_1 * sqrt(1). sqrt is only implemented for c_1 * sqrt(1) and it is not checked that c_1 is positive.
== is implemented, but the other comparison operators are not: casting to a RealIntervalField is implemented instead and the user can compare the intervals.
The == operator is implemented by first reducing A == B to D == 0 and then converting to a different data type (_FactorizedSqrtLinCombination) that can represent linear combinations:
D = c_1 * sqrt(r_{1,1}) * sqrt(r_{1,2}) * ... * sqrt(r_{1,k_1}) + c_2 * sqrt(r_{2,1}) * sqrt(r_{2,2}) * ... * sqrt(r_{2,k_2}) + ... + c_n * sqrt(r_{n,1}) * sqrt(r_{n,2}) * ... * sqrt(r_{n,k_n})
- by just trivially setting
k_i = 0 when r_i = 1 and r_{i,1} = r_i and k_1 = 1 otherwise.
For this data type, multiplying two sqrt(r_{i,j}) with equal r_{i,j} will cancel the two sqrt’s and apply the common r_{i,j} to the c_i of the result instead. Thus, the following procedure for determining whether D == 0 will eventually terminate:
if the number of terms n is 0, return True
if the number of terms n is 1, return c_1 == 0
if there is a r_{i,j} common to each summand, factor it out
pick one of the r_{i,j}, split the sum into two parts “left”, respectively, “right” of all the terms containing sqrt(r_{i,j}), respectively, not containing sqrt(r_{i,j}).
If left^2 - right^2 == 0 is False, return False. (sqrt(r_{i,j})^2 simplifies to r_{i,j} and disappears, so the resulting expression is easier and this recursion terminates eventually.)
If left == 0 (some comment applies), return True
Use interval arithmetic of increasing precision until it is high enough to determine the signs of left and right. Return True if and only if the signs differ, otherwise False.
Examples:
sage: from sage.rings.number_field.number_field import NumberField sage: from sage.rings.integer import Integer sage: from sage.rings.rational import Rational sage: from sage.rings.real_mpfr import RealLiteral, RealField sage: from sage.rings.real_mpfi import RealIntervalField sage: from sage.calculus.var import var sage: from sage.functions.other import sqrt sage: x = var('x') sage: poly = x ** 6 + Rational((3,2))*x**4 + Rational((9,16))*x**2 - Rational((23,64)) sage: nf = NumberField(poly, 'z', embedding = RealField()(0.56227951206)) sage: z = nf.gen() sage: A = SqrtLinCombination(z) sage: B = SqrtLinCombination(Rational((8,9))*z**4 + Rational((10,9))*z**2 + Rational((2,9))) sage: C = SqrtLinCombination(3) sage: D = SqrtLinCombination(Integer(5)) sage: E = SqrtLinCombination(Rational((6,7))) sage: A + B (8/9*z^4 + 10/9*z^2 + z + 2/9) * sqrt(1) sage: B - E (8/9*z^4 + 10/9*z^2 - 40/63) * sqrt(1) sage: A + sqrt(B) * sqrt(B) (8/9*z^4 + 10/9*z^2 + z + 2/9) * sqrt(1) sage: A + sqrt(B) * sqrt(B) + C == A + B + C True sage: A / E (7/6*z) * sqrt(1) sage: B / A.sqrt() (128/207*z^5 + 376/207*z^3 + 302/207*z) * sqrt(z) sage: B / (D * A.sqrt()) (128/1035*z^5 + 376/1035*z^3 + 302/1035*z) * sqrt(z) sage: RIF = RealIntervalField(100) sage: RIF(B.sqrt() + E.sqrt()) 1.73967449622339881238507307209? sage: A - B == 0 False sage: (A + B).sqrt() (1) * sqrt(8/9*z^4 + 10/9*z^2 + z + 2/9) sage: 3 * A.sqrt() + (4 * B).sqrt() + C + 8 == (9 * A).sqrt() + 2 * B.sqrt() + (C * C).sqrt() + 11 - 3 True
- sign()
Returns the +1, 0, -1 depending on whether the value is positive, zero or negative. For the zero case, exact arithmetic is used to certify. Otherwise, interval arithmetic is used.
- sign_with_interval()
Similar to sign, but for the non-zero case, also return the interval certifying the sign - useful for debugging.
Exceptions
All final exceptions are deriving from two base classes:
a subclass of VerifyErrorBase to indicate whether a numerical or exact verification failed
a subclass of EquationType to indicate the type of equation of inequality for which the verification failed.
Intermediate subclasses (those without __init__) are not supposed to be raised.
The hierarchy is as follows:
VerifyErrorBase(RuntimeError)
NumericalVerifyError
InequalityNumericalVerifyError
LogLiftNumericalVerifyError
ExactVerifyError
IsZeroExactVerifyError
EquationType
EdgeEquationType
EdgeEquationExactVerifyError
EdgeEquationLogLiftNumericalVerifyError
CuspConsistencyType
CuspEquationType
CuspEquationExactVerifyError
CuspEquationLogLiftNumericalVerifyError
TiltType
TiltInequalityNumericalVerifyError
TiltProvenPositiveNumericalVerifyError
TiltIsZeroExactVerifyError
ShapeType
ShapePositiveImaginaryPartNumericalVerifyError
- class snappy.verify.exceptions.CuspConsistencyType
A base class indicating that verificatin of an equation involving a cusp failed.
- exception snappy.verify.exceptions.CuspEquationExactVerifyError(value, expected_value)
Exception for failed verification of a polynomial cusp gluing equation using exact arithmetics.
- exception snappy.verify.exceptions.CuspEquationLogLiftNumericalVerifyError(value, expected_value)
Exception for failed numerical verification that a logarithmic cusp equation has error bound by epsilon.
- class snappy.verify.exceptions.CuspEquationType
A base class indicating that a cusp gluing equation (involving the shapes) failed.
- exception snappy.verify.exceptions.EdgeEquationExactVerifyError(value)
Exception for failed verification of a polynomial edge equation using exact arithmetics.
- exception snappy.verify.exceptions.EdgeEquationLogLiftNumericalVerifyError(value)
Exception for failed numerical verification that a logarithmic edge equation has error bound by epsilon.
- class snappy.verify.exceptions.EdgeEquationType
A base class indicating that an edge equation could not be verified.
- class snappy.verify.exceptions.EquationType
A base class to derive subclasses which indicate what kind of equation failed to be verified.
- exception snappy.verify.exceptions.ExactVerifyError
The base for all exceptions resulting from a failed verification of an equation using exact arithmetics.
- exception snappy.verify.exceptions.InequalityNumericalVerifyError
The base for all exceptions resulting from a failed numerical verification of an inequality (typically by interval arithmetics).
- exception snappy.verify.exceptions.IsZeroExactVerifyError
The base for all exceptions resulting from verifying that a desired quantity is zero using exact arithmetics.
- exception snappy.verify.exceptions.LogLiftNumericalVerifyError
To verify a logarithmic gluing equation, the verify module will usually first verify the corresponding polynomial gluing equation. This means that the logarithmic gluing equation will be fulfilled up to a multiple of 2 Pi I. It then computes the logarithms and numerically checks that the result is close (by some epsilon) to the right value. Because we already know that the difference is a multiple of 2 Pi I, checking closeness is enough.
This exception is supposed to be raised if the polynomial gluing equations have passed but checking the logarithmic equation is epsilon-close has failed.
- exception snappy.verify.exceptions.NumericalVerifyError
The base for all exceptions resulting from a failed numerical verification of an equality (using some epsilon) or inequality (typically by interval arithmetics).
- exception snappy.verify.exceptions.ShapePositiveImaginaryPartNumericalVerifyError(value)
Failed numerical verification of a shape having positive imaginary part.
- class snappy.verify.exceptions.ShapeType
Base class for failed verification of legal shapes.
- exception snappy.verify.exceptions.TiltInequalityNumericalVerifyError(value)
Numerically verifying that a tilt is negative has failed.
- exception snappy.verify.exceptions.TiltIsZeroExactVerifyError(value)
Verifying that a tilt is zero has failed using exact arithmetic.
- exception snappy.verify.exceptions.TiltProvenPositiveNumericalVerifyError(value)
Numerically verifying that a tilt is negative has not only failed, we proved that the tilt is positive and thus that this cannot be a proto-canonical triangulation.
- class snappy.verify.exceptions.TiltType
A base class relating to tilts.
- exception snappy.verify.exceptions.VerifyErrorBase
The base for all exceptions related to verification.