Пример #1
0
def _volume_from_shape(z):
    """
    Computes the Bloch-Wigner dilogarithm for z which gives the volume of a
    tetrahedron of the given shape.
    """

    if _within_sage:
        CIF = z.parent()
        if is_ComplexIntervalField(CIF):
            # A different bug in sage:
            # Depending on the sage version, an element in a
            # ComplexIntervalField wouldn't support dilog/polylog, or, even
            # worse, would convert the element to ComplexField first!!!
            #
            # Thus, we convert to ComplexBallField here since the arblib
            # supports a verified interval polylog (albeit giving an interval
            # that seems to be 300 times larger than necessary).
            
            CBF = ComplexBallField(CIF.precision())
            RIF = RealIntervalField(CIF.precision())
    
            return RIF(_unprotected_volume_from_shape(CBF(z)))
        else:
            z = Number(z)

    # Use implementation in number.py that overcomes the cypari bug that you
    # have to explicitly give a precision to dilog, otherwise you lose
    # precision.
    return z.volume()
Пример #2
0
    def _sign_numerical(self, prec):
        """
        Use interval arithmetics with precision prec to try to determine the
        sign. If we could not certify the sign, return None.
        The result is a pair (sign, interval).
        """

        # Evaluate as interval
        RIF = RealIntervalField(prec)
        try:
            interval_val = RIF(self)
        except _SqrtException:
            # This exception happens if we try to take the square root of an
            # interval that contains negative numbers.
            # This is not supposed to happen but if we take the square of a small
            # number and the precision is low, it might happen.
            # It just means we need to use higher precision.
            # So just return "None" to indicate failed certification.
            return None, None

        # Interval certifies positive sign
        if interval_val > 0:
            return +1, interval_val
        # Interval certified negative sign
        if interval_val < 0:
            return -1, interval_val
        # Interval contains zero and we can't say.
        return None, interval_val
Пример #3
0
 def _known_bound(self, iv, post_transform):
     post_transform = normalize_post_transform(self.dop, post_transform)
     Balls = iv.parent()
     Ivs = RealIntervalField(Balls.precision())
     mid = [
         c for c in self._polys.keys()
         if Balls(c).add_error(self._rad(c)).overlaps(iv)
     ]
     mid.sort()
     rad = [self._rad(c) for c in mid]
     crude_bound = (Balls(AnInfinity())
                    if self._is_everywhere_defined() else Balls('nan'))
     if len(mid) < 1 or not mid[0] - rad[0] <= iv <= mid[-1] + rad[-1]:
         return crude_bound
     if not all(
             safe_eq(mid[i] + rad[i], mid[i + 1] - rad[i + 1])
             for i in range(len(mid) - 1)):
         return crude_bound
     bound = None
     for c, r in zip(mid, rad):
         if len(self._polys[c]) <= post_transform.order():
             return crude_bound
         polys = [a.pol for a in self._polys[c]]
         dom = Balls(Ivs(Balls(c).add_error(r)).intersection(Ivs(iv)))
         reduced_dom = dom - c
         img = sum(
             ZZ(j).factorial() * coeff(dom) * polys[j](reduced_dom)
             for j, coeff in enumerate(post_transform))
         bound = img if bound is None else bound.union(img)
     return bound
Пример #4
0
 def subdivide(self, threshold=IR(0.6), factor=IR(0.5)):
     # TODO:
     # - support paths passing very close to singular points
     from sage.rings.real_mpfi import RealIntervalField
     new = [self.vert[0]]
     i = 1
     while i < len(self.vert):
         cur, next = new[-1], self.vert[i]
         rad = cur.dist_to_sing()
         dist_to_next = (next.iv() - cur.iv()).abs()
         if (dist_to_next <= threshold * rad if next.is_ordinary() else
             (cur.value == next.value or cur.is_ordinary()
              and dist_to_next <= threshold * next.dist_to_sing())):
             new.append(next)
             i += 1
         else:
             dir = (next.iv() - cur.iv()) / dist_to_next
             interm = cur.iv() + factor * rad * dir
             is_real = interm.imag().is_zero()
             interm = interm.add_error(rad / 8)
             Step(cur, Point(interm, self.dop)).check_singularity()  # TBI
             my_RIF = RealIntervalField(interm.real().parent().precision())
             if is_real:
                 interm = my_RIF(interm.real()).simplest_rational()
             else:
                 interm = QQi([
                     my_RIF(interm.real()).simplest_rational(),
                     my_RIF(interm.imag()).simplest_rational()
                 ])
             new.append(Point(interm, self.dop))
             logger.debug("subdividing %s -> %s", cur, next)
     new = Path(new, self.dop)
     return new
Пример #5
0
def _test_result(number_field, prec=53, epsilon=1e-10):
    """
        sage: from sage.rings.complex_field import ComplexField
        sage: CF = ComplexField()
        sage: x = var('x')
        sage: nf = NumberField(x**2 + 1, 'x', embedding = CF(1.0j))
        sage: _test_result(nf)
        sage: nf = NumberField(x**2 + 7, 'x', embedding = CF(2.64575j))
        sage: _test_result(nf)
        sage: nf = NumberField(x**8 + 6 * x ** 4 + x + 23, 'x', embedding = CF(0.7747 + 1.25937j))
        sage: _test_result(nf, 212, epsilon = 1e-30)
    """

    CIF = ComplexIntervalField(prec)
    RIF = RealIntervalField(prec)

    real_number_field, x_expression, y_expression = (
        field_containing_real_and_imaginary_part_of_number_field(number_field))

    x_val = x_expression.lift()(RIF(real_number_field.gen_embedding()))
    y_val = y_expression.lift()(RIF(real_number_field.gen_embedding()))
    z_val = CIF(x_val, y_val)

    diff = z_val - CIF(number_field.gen_embedding())

    if not abs(diff) < epsilon:
        raise Exception("Test failed")
Пример #6
0
    def _richcmp_(self, other, op):
        r"""
        Comparison (using the complex embedding).

        TESTS::

            sage: UCF = UniversalCyclotomicField()
            sage: l = [UCF.gen(3), UCF.gen(3)+1, UCF.gen(5), UCF.gen(5,2),
            ....:      UCF.gen(4), 2*UCF.gen(4), UCF.gen(5)-22/3]
            sage: lQQbar = list(map(QQbar,l))
            sage: lQQbar.sort()
            sage: l.sort()
            sage: lQQbar == list(map(QQbar,l))
            True

            sage: for i in range(len(l)):
            ....:     assert l[i] >= l[i] and l[i] <= l[i]
            ....:     for j in range(i):
            ....:         assert l[i] > l[j] and l[j] < l[i]

            sage: fibonacci(200)*(E(5)+E(5,4)) <= fibonacci(199)
            True
            sage: fibonacci(201)*(E(5)+E(5,4)) <= fibonacci(200)
            False
        """
        if self._obj == other._obj:
            return rich_to_bool(op, 0)

        s = self.real_part()
        o = other.real_part()
        if s == o:
            s = self.imag_part()
            o = other.imag_part()

        from sage.rings.real_mpfi import RealIntervalField
        prec = 53
        R = RealIntervalField(prec)
        sa = s._eval_real_(R)
        oa = o._eval_real_(R)
        while sa.overlaps(oa):
            prec <<= 2
            R = RealIntervalField(prec)
            sa = s._eval_real_(R)
            oa = o._eval_real_(R)
        return sa._richcmp_(oa, op)
Пример #7
0
def _rationalize(civ, real=False):
    from sage.rings.real_mpfi import RealIntervalField
    my_RIF = RealIntervalField(civ.real().parent().precision())
    # XXX why can't we do this automatically when civ.imag().contains_zero()???
    if real or civ.imag().is_zero():
        return my_RIF(civ.real()).simplest_rational()
    else:
        return QQi([my_RIF(civ.real()).simplest_rational(),
                    my_RIF(civ.imag()).simplest_rational()])
Пример #8
0
def _rationalize(civ, real=False):
    from sage.rings.real_mpfi import RealIntervalField
    my_RIF = RealIntervalField(civ.real().parent().precision())
    if real or civ.imag().is_zero():
        return my_RIF(civ.real()).simplest_rational()
    else:
        return QQi([
            my_RIF(civ.real()).simplest_rational(),
            my_RIF(civ.imag()).simplest_rational()
        ])
Пример #9
0
def _opposite_signs(left, right, prec):
    """
    Given two objects left and right that can be coerced to real interval of
    the given precision, try to certify their signs. If succeed, return True
    if the signs are opposite and False otherwise. If failed, return None.
    """

    # Try to cast the elements to real intervals
    RIF = RealIntervalField(prec)
    try:
        left_interval = RIF(left)
        right_interval = RIF(right)
    except _SqrtException:
        # This exception happens if we try to take the square root of an
        # interval that contains negative numbers.
        # This is not supposed to happen but if we take the square of a small
        # number and the precision is low, it might happen.
        # It just means we need to use higher precision.
        # So just return "None" to indicate failed certification.
        return None
    
    # Try to determine sign of left expression.
    left_negative    = bool(left_interval  < 0)
    left_positive    = bool(left_interval  > 0)
    left_determined  = left_negative  or left_positive
    
    # Try to determine sign of right expression
    right_negative   = bool(right_interval < 0)
    right_positive   = bool(right_interval > 0)
    right_determined = right_negative or right_positive
    
    # If both signs could be determined
    if left_determined and right_determined:
        # Return true if and only if signs are opposite
        return left_positive ^ right_positive
    
    # At least one sign couldn't be determined.
    return None
Пример #10
0
    def __init__(self, M, initial_shapes, bits_prec=None, dec_prec=None):
        """
        Initializes the KrawczykShapesEngine given an orientable SnapPy
        Manifold M, approximated solutions initial_shapes to the
        gluing equations (e.g., as returned by M.tetrahedra_shapes('rect'))
        and the precision to be used for the desired computations in either
        bits bits_prec or decimal digits dec_prec.

        This requires Sage since it uses Sage's ComplexIntervalField for its
        computations.

        Note that this will choose an independent set of edge equations and
        one equation per cusp. It is known that a solution to such a subset of
        rectangular gluing equations is also a solution to the full set of
        rectangular gluing equations::

            sage: from snappy import Manifold
            sage: M = Manifold("m019")

            sage: C = KrawczykShapesEngine(M, M.tetrahedra_shapes('rect'), bits_prec = 53)
            sage: C.expand_until_certified()
            True
            sage: C.certified_shapes # doctest: +ELLIPSIS
            (0.780552527850...? + 0.914473662967...?*I, 0.780552527850...? + 0.91447366296773?*I, 0.4600211755737...? + 0.6326241936052...?*I)

        Does not work with non-orientable manifolds::

            sage: M = Manifold("m000")
            sage: KrawczykShapesEngine(M, M.tetrahedra_shapes('rect'), bits_prec = 53)
            Traceback (most recent call last):
            ...
            Exception: Manifold needs to be orientable


        Or some non-hyperbolic manifolds::
        
            sage: Manifold("t02333(1,0)").tetrahedra_shapes(intervals = True)
            Traceback (most recent call last):
            ...
            RuntimeError: Could not certify shape intervals, either there are degenerate shapes or the precision must be increased.

        """

        # Require sage
        if not _within_sage:
            raise SageNotAvailable(
                "Sorry, the verify module can only be used within Sage")

        # Convert to precision in bits if necessary
        if dec_prec:
            self.prec = prec_dec_to_bits(dec_prec)
        elif bits_prec:
            self.prec = bits_prec
        else:
            raise Exception("Need dec_prec or bits_prec")

        # Setup interval types of desired precision
        self.CIF = ComplexIntervalField(self.prec)
        self.RIF = RealIntervalField(self.prec)

        # Verify that manifold is orientable
        if not M.is_orientable():
            raise Exception("Manifold needs to be orientable")

        # Initialize the shape intervals, they have zero length
        self.initial_shapes = vector(
            [self.CIF(shape) for shape in initial_shapes])

        # Get an independent set of gluing equations from snap
        self.equations = snap.shapes.enough_gluing_equations(M)
        self._make_sparse_equations()

        self.identity = matrix.identity(self.CIF, len(self.initial_shapes))

        CDF = ComplexDoubleField()

        # Could be sparse
        approx_deriv = self.log_gluing_LHS_derivatives(
            [CDF(shape) for shape in initial_shapes])
        approx_inverse_double = approx_deriv.inverse()
        self.approx_inverse = approx_inverse_double.change_ring(self.CIF)

        # Compute the term z0 - c * f(z0) in the formula for
        # the Krawczyk interval K(z0, [z], f)

        value_at_initial_shapes = self.log_gluing_LHSs(self.initial_shapes)

        self.first_term = (self.initial_shapes -
                           self.approx_inverse * value_at_initial_shapes)

        # Shapes have not been certified yet
        self.certified_shapes = None
Пример #11
0
def __reduce__RealIntervalField(prec, sci_not):
    return RealIntervalField(prec, sci_not)
Пример #12
0
def refine_root(ip, ipd, irt, fld):
    """
    We are given a polynomial and its derivative (with complex
    interval coefficients), an estimated root, and a complex interval
    field to use in computations.  We use interval arithmetic to
    refine the root and prove that we have in fact isolated a unique
    root.

    If we succeed, we return the isolated root; if we fail, we return
    None.

    EXAMPLES::

        sage: from sage.rings.polynomial.complex_roots import *
        sage: x = polygen(ZZ)
        sage: p = x^9 - 1
        sage: ip = CIF['x'](p); ip
        x^9 - 1
        sage: ipd = CIF['x'](p.derivative()); ipd
        9*x^8
        sage: irt = CIF(CC(cos(2*pi/9), sin(2*pi/9))); irt
        0.76604444311897802? + 0.64278760968653926?*I
        sage: ip(irt)
        0.?e-14 + 0.?e-14*I
        sage: ipd(irt)
        6.89439998807080? - 5.78508848717885?*I
        sage: refine_root(ip, ipd, irt, CIF)
        0.766044443118978? + 0.642787609686540?*I
    """

    # There has got to be a better way to do this, but I don't know
    # what it is...

    # We start with a basic fact: if we do an interval Newton-Raphson
    # step, and the refined interval is contained in the original interval,
    # then the refined interval contains exactly one root.

    # Unfortunately, our initial estimated root almost certainly does not
    # contain the actual root (our initial interval is a point, which
    # is exactly equal to whatever floating-point estimate we got from
    # the external solver).  So we need to do multiple Newton-Raphson
    # steps, and check this inclusion property on each step.

    # After a few steps of refinement, if the initial root estimate was
    # close to a root, we should have an essentially perfect interval
    # bound on the root (since Newton-Raphson has quadratic convergence),
    # unless either the real or imaginary component of the root is zero.
    # If the real or imaginary component is zero, then we could spend
    # a long time computing closer and closer approximations to that
    # component.  (This doesn't happen for non-zero components, because
    # of the imprecision of floating-point numbers combined with the
    # outward interval rounding; but close to zero, MPFI provides
    # extremely precise numbers.)

    # If the root is actually a real root, but we start with an imaginary
    # component, we can bounce back and forth between having a positive
    # and negative imaginary component, without ever hitting zero.
    # To deal with this, on every other Newton-Raphson step, instead of
    # replacing the old interval with the new one, we take the union.

    # If the containment check continues to fail many times in a row,
    # we give up and return None; we also return None if we detect
    # that the slope in our current interval is not bounded away
    # from zero at any step.

    # After every refinement step, we check to see if the real or
    # imaginary component of our interval includes zero.  If so, we
    # try setting it to exactly zero.  This gives us a good chance of
    # detecting real roots.  However, we do this replacement at most
    # once per component.

    refinement_steps = 10

    smashed_real = False
    smashed_imag = False

    for i in range(refinement_steps):
        slope = ipd(irt)
        if slope.contains_zero():
            return None
        center = fld(irt.center())
        val = ip(center)

        nirt = center - val / slope
        # print irt, nirt, (nirt in irt), nirt.diameter(), irt.diameter(), center, val, slope
        if nirt in irt and (nirt.diameter() >= irt.diameter() >> 3 or i >= 8):
            # If the new diameter is much less than the original diameter,
            # then we have not yet converged.  (Perhaps we were asked
            # for a particularly high-precision result.)  So we don't
            # return yet.
            return nirt

        if i & 1:
            irt = nirt
        else:
            irt = irt.union(nirt)
            # If we don't find a root after a while, try (approximately)
            # tripling the size of the region.
            if i >= 6:
                rD = irt.real().absolute_diameter()
                iD = irt.imag().absolute_diameter()
                md = max(rD, iD)
                md_intv = RealIntervalField(rD.prec())(-md, md)
                md_cintv = ComplexIntervalField(rD.prec())(md_intv, md_intv)
                irt = irt + md_cintv

        if not smashed_real and irt.real().contains_zero():
            irt = irt.parent()(0, irt.imag())
            smashed_real = True
        if not smashed_imag and irt.imag().contains_zero():
            irt = irt.parent()(irt.real(), 0)
            smashed_imag = True

    return None