Changeset 7425:c26b239c2e4e


Ignore:
Timestamp:
11/25/07 08:10:06 (6 years ago)
Author:
Carl Witty <cwitty@…>
Branch:
default
Message:

Add complex intervals.

Files:
3 added
7 edited

Legend:

Unmodified
Added
Removed
  • sage/rings/all.py

    r6342 r7425  
    103103from complex_number import (is_ComplexNumber, create_ComplexNumber as ComplexNumber) 
    104104Complexes = ComplexField 
     105from complex_interval_field import ComplexIntervalField, is_ComplexIntervalField 
     106from complex_interval import (is_ComplexIntervalFieldElement, create_ComplexIntervalFieldElement as ComplexIntervalFieldElement) 
    105107 
    106108from complex_double import ComplexDoubleField, ComplexDoubleElement, CDF, is_ComplexDoubleElement 
     
    146148 
    147149CC = ComplexField() 
     150CIF = ComplexIntervalField() 
    148151I = CC.gen() 
    149152 
  • sage/rings/complex_field.py

    r4057 r7425  
    2525 
    2626from sage.structure.parent_gens import ParentWithGens 
     27 
     28NumberFieldElement_quadratic = None 
     29def late_import(): 
     30    global NumberFieldElement_quadratic 
     31    if NumberFieldElement_quadratic is None: 
     32        import sage.rings.number_field.number_field_element_quadratic as nfeq 
     33        NumberFieldElement_quadratic = nfeq.NumberFieldElement_quadratic 
    2734 
    2835def is_ComplexField(x): 
     
    7077        0.333333333333333 + 2.00000000000000*I 
    7178 
    72     Note that the second argument is the number of *bits* of precision, 
    73     not the number of digits of precision: 
    74         sage: C(1/3, 2) 
    75         0.333333333333333 + 2.00000000000000*I 
    76  
    7779    We can also coerce rational numbers and integers into C, but 
    78     coercing a polynomial in raising an exception. 
     80    coercing a polynomial will raise an exception. 
    7981 
    8082        sage: Q = RationalField() 
     
    164166            sage: CC(2,3) 
    165167            2.00000000000000 + 3.00000000000000*I 
     168            sage: CC(QQ[I].gen()) 
     169            1.00000000000000*I 
     170            sage: CC.gen() + QQ[I].gen() 
     171            Traceback (most recent call last): 
     172            ... 
     173            TypeError: unsupported operand parent(s) for '+': 'Complex Field with 53 bits of precision' and 'Number Field in I with defining polynomial x^2 + 1' 
    166174        """ 
    167175        if im is None: 
     
    175183                return complex_number.ComplexNumber(self, 
    176184                            sage_eval(x.replace(' ',''), locals={"I":self.gen(),"i":self.gen()})) 
     185 
     186            late_import() 
     187            if isinstance(x, NumberFieldElement_quadratic) and list(x.parent().polynomial()) == [1, 0, 1]: 
     188                (re, im) = list(x) 
     189                return complex_number.ComplexNumber(self, re, im) 
     190 
    177191            try: 
    178192                return x._complex_mpfr_field_( self ) 
     
    266280    def scientific_notation(self, status=None): 
    267281        return self._real_field().scientific_notation(status) 
    268          
    269  
  • sage/rings/complex_number.pyx

    r5809 r7425  
    7373         
    7474        if imag is None: 
     75            if real is None: return 
     76 
    7577            if PY_TYPE_CHECK(real, ComplexNumber): 
    7678                real, imag = (<ComplexNumber>real).real(), (<ComplexNumber>real).imag() 
  • sage/rings/mpfi.pxi

    r4843 r7425  
    213213 
    214214 
    215     int mpfi_has_zero(mpfi_srcptr) 
     215    bint mpfi_has_zero(mpfi_srcptr) 
    216216 
    217217    bint mpfi_nan_p(mpfi_srcptr) 
  • sage/rings/mpfr.pxi

    r6880 r7425  
    114114    #int mpfr_pow_si _PROTO ((mpfr_ptr, mpfr_srcptr, long int, mp_rnd_t));  
    115115 
    116     #int mpfr_min _PROTO ((mpfr_ptr, mpfr_srcptr, mpfr_srcptr, mp_rnd_t)); 
    117     #int mpfr_max _PROTO ((mpfr_ptr, mpfr_srcptr, mpfr_srcptr, mp_rnd_t)); 
     116    int mpfr_min (mpfr_ptr, mpfr_srcptr, mpfr_srcptr, mp_rnd_t) 
     117    int mpfr_max (mpfr_ptr, mpfr_srcptr, mpfr_srcptr, mp_rnd_t) 
    118118 
    119119    int mpfr_fac_ui (mpfr_t rop, unsigned long int op, mp_rnd_t rnd) 
     
    150150    int mpfr_neg (mpfr_ptr rop, mpfr_srcptr op, mp_rnd_t rnd) 
    151151    # int mpfr_eq (mpfr_srcptr rop, mpfr_srcptr op, unsigned long i) 
     152    int mpfr_cmp (mpfr_t op1, mpfr_t op2) 
     153    bint mpfr_greater_p (mpfr_t op1, mpfr_t op2) 
     154    bint mpfr_greaterequal_p (mpfr_t op1, mpfr_t op2) 
    152155    bint mpfr_less_p (mpfr_t op1, mpfr_t op2) 
    153156    bint mpfr_lessequal_p (mpfr_t op1, mpfr_t op2) 
    154     int mpfr_cmp (mpfr_t op1, mpfr_t op2) 
     157    bint mpfr_lessgreater_p (mpfr_t op1, mpfr_t op2) 
    155158    bint mpfr_equal_p (mpfr_t op1, mpfr_t op2) 
     159    bint mpfr_unordered_p (mpfr_t op1, mpfr_t op2) 
  • sage/rings/real_mpfi.pyx

    r7244 r7425  
    106106############################################################################ 
    107107  
    108 #***************************************************************************** 
    109 # general TODOs:  
    110 #  
    111 # more type conversions and coercion. examples: 
    112 # sage: R(1/2) 
    113 # TypeError: Unable to convert x (='1/2') to mpfi. 
    114 #  
    115 # sage: 1 + R(42) 
    116 # _49 = 1 
    117 #***************************************************************************** 
    118  
    119108import math # for log 
    120109import sys 
     
    430419            else: 
    431420                raise TypeError, "Canonical coercion from lower to higher precision not defined" 
    432         if isinstance(x, (Integer, Rational, sage.rings.algebraic_real.AlgebraicRealNumber)): 
     421        if isinstance(x, (Integer, Rational)): 
    433422            return self(x) 
    434423        cdef RealNumber lower, upper 
     
    475464            raise IndexError 
    476465 
    477 #    def complex_field(self): 
    478 #        """ 
    479 #        Return complex field of the same precision. 
    480 #        """ 
    481 #        return sage.rings.complex_field.ComplexField(self.prec()) 
     466    def complex_field(self): 
     467        """ 
     468        Return complex field of the same precision. 
     469 
     470        EXAMPLES: 
     471            sage: RIF.complex_field() 
     472            Complex Interval Field with 53 bits of precision 
     473        """ 
     474        return sage.rings.complex_interval_field.ComplexIntervalField(self.prec()) 
    482475 
    483476    def ngens(self): 
     
    538531        return self.__prec 
    539532     
    540     # int mpfi_const_pi (mpfi_ptr) 
    541533    def pi(self): 
    542534        """ 
     
    558550        return x 
    559551 
    560  
    561     # int mpfi_const_euler (mpfi_ptr) 
    562552    def euler_constant(self): 
    563553        """ 
     
    10521042        return x 
    10531043 
     1044    def is_exact(self): 
     1045        return mpfr_equal_p(&self.value.left, &self.value.right)             
     1046 
    10541047    def magnitude(self): 
    10551048        """ 
     
    18221815            return mpfr_lessequal_p(&rt.value.right, &lt.value.left) 
    18231816 
     1817    def __nonzero__(self): 
     1818        """ 
     1819        Returns true if self is not known to be exactly zero. 
     1820 
     1821        EXAMPLES: 
     1822            sage: RIF(0).__nonzero__() 
     1823            False 
     1824            sage: RIF(1).__nonzero__() 
     1825            True 
     1826            sage: RIF(1, 2).__nonzero__() 
     1827            True 
     1828            sage: RIF(0, 1).__nonzero__() 
     1829            True 
     1830            sage: RIF(-1, 1).__nonzero__() 
     1831            True 
     1832        """ 
     1833        return not (mpfr_zero_p(&self.value.left) and mpfr_zero_p(&self.value.right)) 
     1834 
    18241835    def __cmp__(left, right): 
    18251836        """ 
     
    19061917        except TypeError, msg: 
    19071918            return False 
     1919 
     1920    def contains_zero(self): 
     1921        """ 
     1922        Returns True if self is an interval containing zero. 
     1923 
     1924        EXAMPLES: 
     1925            sage: RIF(0).contains_zero() 
     1926            True 
     1927            sage: RIF(1, 2).contains_zero() 
     1928            False 
     1929            sage: RIF(-1, 1).contains_zero() 
     1930            True 
     1931            sage: RIF(-1, 0).contains_zero() 
     1932            True 
     1933        """ 
     1934        return mpfi_has_zero(self.value) 
     1935 
     1936    def overlaps(self, RealIntervalFieldElement other): 
     1937        """ 
     1938        Returns true if self and other are intervals with at least 
     1939        one value in common.  For intervals a and b, we have 
     1940        a.overlaps(b) iff not(a!=b). 
     1941 
     1942        EXAMPLES: 
     1943            sage: RIF(0, 1).overlaps(RIF(1, 2)) 
     1944            True 
     1945            sage: RIF(1, 2).overlaps(RIF(0, 1)) 
     1946            True 
     1947            sage: RIF(0, 1).overlaps(RIF(2, 3)) 
     1948            False 
     1949            sage: RIF(2, 3).overlaps(RIF(0, 1)) 
     1950            False 
     1951            sage: RIF(0, 3).overlaps(RIF(1, 2)) 
     1952            True 
     1953            sage: RIF(0, 2).overlaps(RIF(1, 3)) 
     1954            True 
     1955        """ 
     1956        return mpfr_greaterequal_p(&self.value.right, &other.value.left) \ 
     1957           and mpfr_greaterequal_p(&other.value.right, &self.value.left) 
    19081958 
    19091959    def intersection(self, other): 
     
    20962146#             exponent = x._parent(exponent) 
    20972147#         return self.__pow(exponent) 
     2148 
     2149    def __pow__(self, exponent, modulus): 
     2150        if isinstance(exponent, (int, long, Integer)): 
     2151            return sage.rings.ring_element.RingElement.__pow__(self, exponent) 
     2152        return (self.log() * exponent).exp() 
     2153         
    20982154 
    20992155    def log(self, base='e'): 
  • setup.py

    r7388 r7425  
    659659    Extension('sage.rings.real_mpfi', 
    660660              sources = ['sage/rings/real_mpfi.pyx'], 
     661              libraries = ['mpfi', 'mpfr', 'gmp']), \ 
     662 
     663    Extension('sage.rings.complex_interval', 
     664              sources = ['sage/rings/complex_interval.pyx'], 
    661665              libraries = ['mpfi', 'mpfr', 'gmp']), \ 
    662666 
Note: See TracChangeset for help on using the changeset viewer.