Changeset 7425:c26b239c2e4e
- Timestamp:
- 11/25/07 08:10:06 (6 years ago)
- Branch:
- default
- Files:
-
- 3 added
- 7 edited
-
sage/rings/all.py (modified) (2 diffs)
-
sage/rings/complex_field.py (modified) (5 diffs)
-
sage/rings/complex_interval.pxd (added)
-
sage/rings/complex_interval.pyx (added)
-
sage/rings/complex_interval_field.py (added)
-
sage/rings/complex_number.pyx (modified) (1 diff)
-
sage/rings/mpfi.pxi (modified) (1 diff)
-
sage/rings/mpfr.pxi (modified) (2 diffs)
-
sage/rings/real_mpfi.pyx (modified) (9 diffs)
-
setup.py (modified) (1 diff)
Legend:
- Unmodified
- Added
- Removed
-
sage/rings/all.py
r6342 r7425 103 103 from complex_number import (is_ComplexNumber, create_ComplexNumber as ComplexNumber) 104 104 Complexes = ComplexField 105 from complex_interval_field import ComplexIntervalField, is_ComplexIntervalField 106 from complex_interval import (is_ComplexIntervalFieldElement, create_ComplexIntervalFieldElement as ComplexIntervalFieldElement) 105 107 106 108 from complex_double import ComplexDoubleField, ComplexDoubleElement, CDF, is_ComplexDoubleElement … … 146 148 147 149 CC = ComplexField() 150 CIF = ComplexIntervalField() 148 151 I = CC.gen() 149 152 -
sage/rings/complex_field.py
r4057 r7425 25 25 26 26 from sage.structure.parent_gens import ParentWithGens 27 28 NumberFieldElement_quadratic = None 29 def 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 27 34 28 35 def is_ComplexField(x): … … 70 77 0.333333333333333 + 2.00000000000000*I 71 78 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*I76 77 79 We can also coerce rational numbers and integers into C, but 78 coercing a polynomial in raisingan exception.80 coercing a polynomial will raise an exception. 79 81 80 82 sage: Q = RationalField() … … 164 166 sage: CC(2,3) 165 167 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' 166 174 """ 167 175 if im is None: … … 175 183 return complex_number.ComplexNumber(self, 176 184 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 177 191 try: 178 192 return x._complex_mpfr_field_( self ) … … 266 280 def scientific_notation(self, status=None): 267 281 return self._real_field().scientific_notation(status) 268 269 -
sage/rings/complex_number.pyx
r5809 r7425 73 73 74 74 if imag is None: 75 if real is None: return 76 75 77 if PY_TYPE_CHECK(real, ComplexNumber): 76 78 real, imag = (<ComplexNumber>real).real(), (<ComplexNumber>real).imag() -
sage/rings/mpfi.pxi
r4843 r7425 213 213 214 214 215 int mpfi_has_zero(mpfi_srcptr)215 bint mpfi_has_zero(mpfi_srcptr) 216 216 217 217 bint mpfi_nan_p(mpfi_srcptr) -
sage/rings/mpfr.pxi
r6880 r7425 114 114 #int mpfr_pow_si _PROTO ((mpfr_ptr, mpfr_srcptr, long int, mp_rnd_t)); 115 115 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) 118 118 119 119 int mpfr_fac_ui (mpfr_t rop, unsigned long int op, mp_rnd_t rnd) … … 150 150 int mpfr_neg (mpfr_ptr rop, mpfr_srcptr op, mp_rnd_t rnd) 151 151 # 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) 152 155 bint mpfr_less_p (mpfr_t op1, mpfr_t op2) 153 156 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) 155 158 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 106 106 ############################################################################ 107 107 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 = 1117 #*****************************************************************************118 119 108 import math # for log 120 109 import sys … … 430 419 else: 431 420 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)): 433 422 return self(x) 434 423 cdef RealNumber lower, upper … … 475 464 raise IndexError 476 465 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()) 482 475 483 476 def ngens(self): … … 538 531 return self.__prec 539 532 540 # int mpfi_const_pi (mpfi_ptr)541 533 def pi(self): 542 534 """ … … 558 550 return x 559 551 560 561 # int mpfi_const_euler (mpfi_ptr)562 552 def euler_constant(self): 563 553 """ … … 1052 1042 return x 1053 1043 1044 def is_exact(self): 1045 return mpfr_equal_p(&self.value.left, &self.value.right) 1046 1054 1047 def magnitude(self): 1055 1048 """ … … 1822 1815 return mpfr_lessequal_p(&rt.value.right, <.value.left) 1823 1816 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 1824 1835 def __cmp__(left, right): 1825 1836 """ … … 1906 1917 except TypeError, msg: 1907 1918 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) 1908 1958 1909 1959 def intersection(self, other): … … 2096 2146 # exponent = x._parent(exponent) 2097 2147 # 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 2098 2154 2099 2155 def log(self, base='e'): -
setup.py
r7388 r7425 659 659 Extension('sage.rings.real_mpfi', 660 660 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'], 661 665 libraries = ['mpfi', 'mpfr', 'gmp']), \ 662 666
Note: See TracChangeset
for help on using the changeset viewer.
