Skip to content

Commit

Permalink
Trac #19330: Implement conversion of interval fields to real/complex …
Browse files Browse the repository at this point in the history
…fields

It surprises me that this isn't implemented yet:
{{{
sage: RR(RIF("1.2"))
TypeError: Unable to convert x (='1.200000000000000?') to real number.
}}}

URL: http://trac.sagemath.org/19330
Reported by: jdemeyer
Ticket author(s): Jeroen Demeyer
Reviewer(s): Vincent Delecroix
  • Loading branch information
Release Manager authored and vbraun committed Oct 12, 2015
2 parents a463174 + 752c401 commit 2cb030d
Show file tree
Hide file tree
Showing 5 changed files with 87 additions and 69 deletions.
22 changes: 19 additions & 3 deletions src/sage/rings/complex_interval.pyx
Original file line number Diff line number Diff line change
Expand Up @@ -51,10 +51,9 @@ from complex_number cimport ComplexNumber
import complex_interval_field
from complex_field import ComplexField
import sage.misc.misc
import integer
cimport integer
import infinity
import real_mpfi
import real_mpfr
cimport real_mpfi
cimport real_mpfr


Expand Down Expand Up @@ -927,6 +926,23 @@ cdef class ComplexIntervalFieldElement(sage.structure.element.FieldElement):

return x

def _complex_mpfr_field_(self, field):
"""
Convert to a complex field.
EXAMPLES::
sage: re = RIF("1.2")
sage: im = RIF(2, 3)
sage: a = ComplexIntervalField(30)(re, im)
sage: CC(a)
1.20000000018626 + 2.50000000000000*I
"""
cdef ComplexNumber x = field(0)
mpfi_mid(x.__re, self.__re)
mpfi_mid(x.__im, self.__im)
return x

def __int__(self):
"""
Convert ``self`` to an ``int``.
Expand Down
9 changes: 7 additions & 2 deletions src/sage/rings/complex_interval_field.py
Original file line number Diff line number Diff line change
Expand Up @@ -377,10 +377,15 @@ def __call__(self, x, im=None):
2 + 3*I
sage: CIF(pi, e)
3.141592653589794? + 2.718281828459046?*I
sage: ComplexIntervalField(100)(CIF(RIF(2,3)))
3.?
"""
if im is None:
if isinstance(x, complex_interval.ComplexIntervalFieldElement) and x.parent() is self:
return x
if isinstance(x, complex_interval.ComplexIntervalFieldElement):
if x.parent() is self:
return x
else:
return complex_interval.ComplexIntervalFieldElement(self, x)
elif isinstance(x, complex_double.ComplexDoubleElement):
return complex_interval.ComplexIntervalFieldElement(self, x.real(), x.imag())
elif isinstance(x, str):
Expand Down
26 changes: 6 additions & 20 deletions src/sage/rings/qqbar.py
Original file line number Diff line number Diff line change
Expand Up @@ -362,8 +362,8 @@
AA(2)
Just for fun, let's try ``sage_input`` on a very complicated expression. The
output of this example changed with the rewritting of polynomial multiplication
algorithms in #10255::
output of this example changed with the rewriting of polynomial multiplication
algorithms in :trac:`10255`::
sage: rt2 = sqrt(AA(2))
sage: rt3 = sqrt(QQbar(3))
Expand Down Expand Up @@ -4553,13 +4553,13 @@ def complex_number(self, field):
EXAMPLES::
sage: a = QQbar.zeta(5)
sage: a.complex_number(CIF)
sage: a.complex_number(CC)
0.309016994374947 + 0.951056516295154*I
sage: (a + a.conjugate()).complex_number(CIF)
sage: (a + a.conjugate()).complex_number(CC)
0.618033988749895 - 5.42101086242752e-20*I
"""
v = self.interval(ComplexIntervalField(field.prec()))
return v.center()
return field(v)

def complex_exact(self, field):
r"""
Expand Down Expand Up @@ -5204,21 +5204,7 @@ def real_number(self, field):
1.41421356237309
"""
v = self.interval(RealIntervalField(field.prec()))

mode = field.rounding_mode()
if mode == 'RNDN':
return v.center()
if mode == 'RNDD':
return v.lower()
if mode == 'RNDU':
return v.upper()
if mode == 'RNDZ':
if v > 0:
return field(v.lower())
elif v < 0:
return field(v.upper())
else:
return field(0)
return field(v)

_mpfr_ = real_number

Expand Down
6 changes: 2 additions & 4 deletions src/sage/rings/real_mpfi.pxd
Original file line number Diff line number Diff line change
Expand Up @@ -2,15 +2,13 @@ from sage.libs.mpfi cimport *

cimport sage.rings.ring

cimport sage.structure.element
from sage.structure.element cimport RingElement

from rational import Rational
from rational cimport Rational

cimport real_mpfr

cdef class RealIntervalFieldElement(sage.structure.element.RingElement) # forward decl
cdef class RealIntervalFieldElement(RingElement) # forward decl

cdef class RealIntervalField_class(sage.rings.ring.Field):
cdef int __prec
Expand All @@ -35,7 +33,7 @@ cdef class RealIntervalField_class(sage.rings.ring.Field):
cdef RealIntervalFieldElement _new(self)


cdef class RealIntervalFieldElement(sage.structure.element.RingElement):
cdef class RealIntervalFieldElement(RingElement):
cdef mpfi_t value
cdef char init
cdef RealIntervalFieldElement _new(self)
Expand Down
93 changes: 53 additions & 40 deletions src/sage/rings/real_mpfi.pyx
Original file line number Diff line number Diff line change
Expand Up @@ -239,39 +239,29 @@ Comparisons with numpy types are right (see :trac:`17758` and :trac:`18076`)::

import math # for log
import sys
import operator

include 'sage/ext/interrupt.pxi'
include "sage/ext/cdefs.pxi"
from cpython.mem cimport *
from cpython.string cimport *

cimport sage.rings.ring
import sage.rings.ring

cimport sage.structure.element
from sage.structure.element cimport RingElement, Element, ModuleElement
import sage.structure.element

cimport real_mpfr
from real_mpfr cimport RealField_class, RealNumber
from real_mpfr import RealField
import real_mpfr

import operator
from real_mpfr cimport RealField_class, RealNumber, RealField
from sage.libs.mpfr cimport MPFR_RNDN, MPFR_RNDZ, MPFR_RNDU, MPFR_RNDD, MPFR_RNDA

from integer import Integer
from integer cimport Integer

from real_double import RealDoubleElement
from real_double cimport RealDoubleElement

import sage.rings.complex_field

import sage.rings.infinity

from sage.structure.parent_gens cimport ParentWithGens

cdef class RealIntervalFieldElement(sage.structure.element.RingElement)

#*****************************************************************************
#
Expand Down Expand Up @@ -1115,16 +1105,13 @@ cdef class RealIntervalField_class(sage.rings.ring.Field):
return self(-1)
raise ValueError, "No %sth root of unity in self"%n

R = RealIntervalField()

#*****************************************************************************
#
# RealIntervalFieldElement -- element of Real Field
#
#
#
#*****************************************************************************
cdef class RealIntervalFieldElement(sage.structure.element.RingElement):
cdef class RealIntervalFieldElement(RingElement):
"""
A real number interval.
"""
Expand Down Expand Up @@ -3073,33 +3060,59 @@ cdef class RealIntervalFieldElement(sage.structure.element.RingElement):
# Conversions
###########################################

# def __float__(self):
# return mpfr_get_d(self.value, (<RealIntervalField>self._parent).rnd)

# def __int__(self):
# """
# Returns integer truncation of this real number.
# """
# s = self.str(32)
# i = s.find('.')
# return int(s[:i], 32)
def _mpfr_(self, RealField_class field):
"""
Convert to a real field, honoring the rounding mode of the
real field.
# def __long__(self):
# """
# Returns long integer truncation of this real number.
# """
# s = self.str(32)
# i = s.find('.')
# return long(s[:i], 32)
EXAMPLES::
# def __complex__(self):
# return complex(float(self))
sage: a = RealIntervalField(30)("1.2")
sage: RR(a)
1.20000000018626
sage: b = RIF(-1, 3)
sage: RR(b)
1.00000000000000
# def _complex_number_(self):
# return sage.rings.complex_field.ComplexField(self.prec())(self)
With different rounding modes::
# def _pari_(self):
# return sage.libs.pari.all.pari.new_with_bits_prec(str(self), (<RealIntervalField>self._parent).__prec)
sage: RealField(53, rnd="RNDU")(a)
1.20000000111759
sage: RealField(53, rnd="RNDD")(a)
1.19999999925494
sage: RealField(53, rnd="RNDZ")(a)
1.19999999925494
sage: RealField(53, rnd="RNDU")(b)
3.00000000000000
sage: RealField(53, rnd="RNDD")(b)
-1.00000000000000
sage: RealField(53, rnd="RNDZ")(b)
0.000000000000000
"""
cdef RealNumber x = field._new()
if field.rnd == MPFR_RNDN:
mpfi_mid(x.value, self.value)
elif field.rnd == MPFR_RNDD:
mpfi_get_left(x.value, self.value)
elif field.rnd == MPFR_RNDU:
mpfi_get_right(x.value, self.value)
elif field.rnd == MPFR_RNDZ:
if mpfi_is_strictly_pos_default(self.value): # interval is > 0
mpfi_get_left(x.value, self.value)
elif mpfi_is_strictly_neg_default(self.value): # interval is < 0
mpfi_get_right(x.value, self.value)
else:
mpfr_set_zero(x.value, 1) # interval contains 0
elif field.rnd == MPFR_RNDA:
# return the endpoint which is furthest from 0
lo, hi = self.endpoints()
if hi.abs() >= lo.abs():
mpfi_get_right(x.value, self.value)
else:
mpfi_get_left(x.value, self.value)
else:
raise AssertionError("%s has unknown rounding mode"%field)
return x

def unique_sign(self):
r"""
Expand Down

0 comments on commit 2cb030d

Please sign in to comment.