source: sage/rings/real_mpfr.pyx @ 3284:6d5a1bc98001

Revision 3284:6d5a1bc98001, 60.5 KB checked in by William Stein <wstein@…>, 6 years ago (diff)

merge

Line 
1"""
2Field of Arbitrary Precision Real Numbers
3
4AUTHORS:
5   -- Kyle Schalm <kschalm@math.utexas.edu> (2005-09)
6   -- William Stein <wstein@gmail.com>: bug fixes, examples, maintenance
7   -- Didier Deshommes <dfdeshom@gmail.com> (2006-03-19): examples
8   -- David Harvey (2006-09-20): compatibility with Element._parent
9   -- William Stein (2006-10): default printing truncates to avoid base-2
10              rounding confusing (fix suggested by Bill Hart)
11
12EXAMPLES:
13
14
15
16A difficult conversion:
17
18    sage: RR(sys.maxint)
19    9223372036854770000     # 64-bit
20    2147483647.00000        # 32-bit
21
22TESTS:
23    sage: -1e30
24    -1000000000000000000000000000000
25"""
26 
27#*****************************************************************************
28#
29#   SAGE: System for Algebra and Geometry Experimentation   
30#
31#       Copyright (C) 2005-2006 William Stein <wstein@gmail.com>
32#
33#  Distributed under the terms of the GNU General Public License (GPL)
34#
35#    This code is distributed in the hope that it will be useful,
36#    but WITHOUT ANY WARRANTY; without even the implied warranty of
37#    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
38#    General Public License for more details.
39#
40#  The full text of the GPL is available at:
41#
42#                  http://www.gnu.org/licenses/
43#*****************************************************************************
44 
45import math # for log
46import sys
47
48include '../ext/interrupt.pxi'
49include "../ext/stdsage.pxi"
50
51cimport sage.rings.ring
52import  sage.rings.ring
53
54cimport sage.structure.element
55from sage.structure.element cimport RingElement, Element, ModuleElement
56import  sage.structure.element
57
58import sage.misc.misc as misc
59
60import sage.structure.coerce
61import operator
62
63from sage.libs.pari.gen import PariInstance, gen
64from sage.libs.pari.gen cimport PariInstance, gen
65
66from integer import Integer
67from integer cimport Integer
68from rational import Rational
69from rational cimport Rational
70
71from real_double import is_RealDoubleElement
72
73import sage.rings.complex_field
74
75import sage.rings.infinity
76
77from sage.structure.parent_gens cimport ParentWithGens
78
79#*****************************************************************************
80# Headers.  When you past things in here from mpfr, be sure
81# to remove const's, since those aren't allowed in pyrex.  Also, it can be
82# challenging figuring out how to modify things from mpfr.h to be valid pyrex
83# code.    Note that what is here is only used for generating the C code.
84# The C compiler doesn't see any of this -- it only sees mpfr.h and stdlib.h
85#*****************************************************************************
86
87cdef class RealNumber(sage.structure.element.RingElement)
88
89#*****************************************************************************
90#
91#       Implementation
92#
93#*****************************************************************************
94
95#*****************************************************************************
96#
97#       External Python access to constants
98#
99#*****************************************************************************
100
101def mpfr_prec_min():
102    """
103    Return the mpfr variable MPFR_PREC_MIN.
104    """
105    return MPFR_PREC_MIN
106
107def mpfr_prec_max():
108    return MPFR_PREC_MAX
109
110#*****************************************************************************
111#
112#       Real Field
113#
114#*****************************************************************************
115# The real field is in Pyrex, so mpfr elements will have access to
116# their parent via direct C calls, which will be faster.
117
118_rounding_modes = ['RNDN', 'RNDZ', 'RNDU', 'RNDD']
119
120cdef class RealField(sage.rings.ring.Field):
121    """
122    RealField(prec, sci_not, rnd):
123
124    INPUT:
125        prec -- (integer) precision; default = 53
126                prec is the number of bits used to represent the
127                mantissa of a floating-point number.  The
128                precision can be any integer between mpfr_prec_min()
129                and mpfr_prec_max(). In the current implementation,
130                mpfr_prec_min() is equal to 2.
131
132        sci_not -- (default: False) if True, always display
133                using scientific notation; if False, display
134                using scientific notation only for very large or
135                very small numbers
136
137        rnd -- (string) the rounding mode
138                RNDN -- (default) round to nearest: Knuth says this is
139                        the best choice to prevent ``floating point
140                        drift''.
141                RNDD -- round towards minus infinity
142                RNDZ -- round towards zero
143                RNDU -- round towards plus infinity
144
145    EXAMPLES:
146        sage: RealField(10)
147        Real Field with 10 bits of precision
148        sage: RealField()
149        Real Field with 53 bits of precision
150        sage: RealField(100000)
151        Real Field with 100000 bits of precision
152
153    NOTE: The default precision is 53, since according to the MPFR
154       manual: 'mpfr should be able to exactly reproduce all
155       computations with double-precision machine floating-point
156       numbers (double type in C), except the default exponent
157       range is much wider and subnormal numbers are not
158       implemented.'
159    """
160
161    def __init__(self, int prec=53, int sci_not=0, rnd="RNDN"):
162        if prec < MPFR_PREC_MIN or prec > MPFR_PREC_MAX:
163            raise ValueError, "prec (=%s) must be >= %s and <= %s."%(
164                prec, MPFR_PREC_MIN, MPFR_PREC_MAX)
165        self.__prec = prec
166        if not isinstance(rnd, str):
167            raise TypeError, "rnd must be a string"
168        self.sci_not = sci_not
169        n = _rounding_modes.index(rnd)
170        if n == -1:
171            raise ValueError, "rnd (=%s) must be one of RNDN, RNDZ, RNDU, or RNDD"%rnd
172        self.rnd = n
173        self.rnd_str = rnd
174        ParentWithGens.__init__(self, self, tuple([]), False)
175
176    cdef RealNumber _new(self):
177        """
178        Return a new real number with parent self.
179        """
180        cdef RealNumber x
181        x = PY_NEW(RealNumber)
182        x._parent = self
183        mpfr_init2(x.value, self.__prec)
184        x.init = 1
185        return x
186
187    def _repr_(self):
188        s = "Real Field with %s bits of precision"%self.__prec
189        if self.rnd != GMP_RNDN:
190            s = s + " and rounding %s"%(self.rnd_str)
191        return s
192
193    def _latex_(self):
194        return "\\R"
195
196    def is_exact(self):
197        return False
198
199    def __call__(self, x, base=10):
200        """
201        Coerce x into this real field.
202
203        EXAMPLES:
204            sage: R = RealField(20)
205            sage: R('1.234')
206            1.2339
207            sage: R('2', base=2)
208            Traceback (most recent call last):
209            ...
210            TypeError: Unable to convert x (='2') to real number.
211            sage: a = R('1.1001', base=2); a
212            1.5625
213            sage: a.str(2)
214            '1.1001000000000000000'
215        """
216        if hasattr(x, '_mpfr_'):
217            # This design with the hasattr is very annoying.
218            # The only thing that uses it right now is symbolic constants
219            # and symbolic function evaluation.
220            # Getting rid of this would speed things up.
221            return x._mpfr_(self)
222        cdef RealNumber z
223        z = self._new()
224        z._set(x, base)
225        return z
226
227    cdef _coerce_c_impl(self, x):
228        """
229        Canonical coercion of x to this mpfr real field.
230
231        The rings that canonically coerce to this mpfr real field are:
232             * this real field itself
233             * any other mpfr real field with precision that is as large as this one
234             * int, long, integer, and rational rings.
235             * real mathematical constants
236        """
237        if isinstance(x, RealNumber):
238            P = x.parent()
239            if (<RealField> P).__prec >= self.__prec:
240                return self(x)
241            else:
242                raise TypeError, "Canonical coercion from lower to higher precision not defined"
243        elif isinstance(x, (Integer, Rational)):
244            return self(x)
245        elif self.__prec <= 53 and is_RealDoubleElement(x):
246            return self(x)
247        import sage.functions.constants
248        return self._coerce_try(x, [sage.functions.constants.ConstantRing])
249
250    def __cmp__(self, other):
251        """
252        EXAMPLES:
253            sage: RealField(10) == RealField(11)
254            False
255            sage: RealField(10) == RealField(10)
256            True
257            sage: RealField(10,rnd='RNDN') == RealField(10,rnd='RNDZ')
258            False
259            sage: RealField(10,sci_not=True) == RealField(10,sci_not=False)
260            False
261            sage: RealField(10) == IntegerRing()
262            False
263        """
264        if not isinstance(other, RealField):
265            return -1
266        cdef RealField _other
267        _other = other  # to access C structure
268        if self.__prec == _other.__prec and self.rnd == _other.rnd \
269               and self.sci_not == _other.sci_not:
270            return 0
271        return 1
272
273    def __reduce__(self):
274        """
275        EXAMPLES:
276            sage: R = RealField(sci_not=1, prec=200, rnd='RNDU')
277            sage: loads(dumps(R)) == R
278            True
279        """
280        return __create__RealField_version0, (self.__prec, self.sci_not, self.rnd_str)
281
282    def gen(self, i=0):
283        if i == 0:
284            return self(1)
285        else:
286            raise IndexError
287
288    def complex_field(self):
289        """
290        Return complex field of the same precision.
291        """
292        return sage.rings.complex_field.ComplexField(self.prec())
293
294    def ngens(self):
295        return 1
296
297    def gens(self):
298        return [self.gen()]
299
300    def _is_valid_homomorphism_(self, codomain, im_gens):
301        try:
302            s = codomain._coerce_(self(1))
303        except TypeError:
304            return False
305        return s == im_gens[0]
306
307    def is_atomic_repr(self):
308        """
309        Returns True, to signify that elements of this field print
310        without sums, so parenthesis aren't required, e.g., in
311        coefficients of polynomials.
312
313        EXAMPLES:
314            sage: RealField(10).is_atomic_repr()
315            True
316        """
317        return True
318
319    def is_finite(self):
320        """
321        Returns False, since the field of real numbers is not finite.
322       
323        EXAMPLES:
324            sage: RealField(10).is_finite()
325            False
326        """
327        return False
328
329    def characteristic(self):
330        """
331        Returns 0, since the field of real numbers has characteristic 0.
332
333        EXAMPLES:
334            sage: RealField(10).characteristic()
335            0
336        """
337        return 0
338   
339    def name(self):
340        return "RealField%s_%s"%(self.__prec,self.rnd)
341   
342    def __hash__(self):
343        return hash(self.name())
344
345    def precision(self):
346        return self.__prec
347
348    def prec(self):
349        return self.__prec
350   
351    # int mpfr_const_pi (mpfr_t rop, mp_rnd_t rnd)
352    def pi(self):
353        """
354        Returns pi to the precision of this field.
355
356        EXAMPLES:
357            sage: R = RealField(100)
358            sage: R.pi()
359            3.1415926535897932384626433832
360            sage: R.pi().sqrt()/2
361            0.88622692545275801364908374167
362            sage: R = RealField(150)
363            sage: R.pi().sqrt()/2
364            0.88622692545275801364908374167057259139877472
365        """
366        cdef RealNumber x
367        x = self._new()       
368        mpfr_const_pi(x.value, self.rnd)
369        return x
370
371
372    # int mpfr_const_euler (mpfr_t rop, mp_rnd_t rnd)
373    def euler_constant(self):
374        """
375        Returns Euler's gamma constant to the precision of this field.
376
377        EXAMPLES:
378            sage: RealField(100).euler_constant()
379            0.57721566490153286060651209008
380        """
381        cdef RealNumber x
382        x = self._new()
383        mpfr_const_euler(x.value, self.rnd)
384        return x
385
386    # int mpfr_const_catalan (mpfr_t rop, mp_rnd_t rnd)
387    def catalan_constant(self):
388        """
389        Returns Catalan's constant to the precision of this field.
390
391        EXAMPLES:
392            sage: RealField(100).catalan_constant()
393            0.91596559417721901505460351493
394        """
395        cdef RealNumber x
396        x = self._new()
397        mpfr_const_catalan(x.value, self.rnd)
398        return x
399
400    # int mpfr_const_log2 (mpfr_t rop, mp_rnd_t rnd)
401    def log2(self):
402        """
403        Returns log(2) to the precision of this field.
404
405        EXAMPLES:
406            sage: R=RealField(100)
407            sage: R.log2()
408            0.69314718055994530941723212145
409            sage: R(2).log()
410            0.69314718055994530941723212145
411        """
412        cdef RealNumber x
413        x = self._new()       
414        mpfr_const_log2(x.value, self.rnd)
415        return x
416
417    def factorial(self, int n):
418        """
419        Return the factorial of the integer n as a real number.
420        """
421        cdef RealNumber x
422        if n < 0:
423            raise ArithmeticError, "n must be nonnegative"
424        x = self._new()
425        mpfr_fac_ui(x.value, n, self.rnd)
426        return x
427
428    def rounding_mode(self):
429        return _rounding_modes[self.rnd]
430
431    def scientific_notation(self, status=None):
432        """
433        Set or return the scientific notation printing flag.  If this flag
434        is True then real numbers with this space as parent print using
435        scientific notation.
436       
437        INPUT:
438            status -- (bool --) optional flag
439        """
440        if status is None:
441            return bool(self.sci_not)
442        else:
443            self.sci_not = status
444       
445    def zeta(self, n=2):
446        """
447        Return an $n$-th root of unity in the real field,
448        if one exists, or raise a ValueError otherwise.
449
450        EXAMPLES:
451            sage: R = RealField()
452            sage: R.zeta()
453            -1.00000000000000
454            sage: R.zeta(1)
455            1.00000000000000
456            sage: R.zeta(5)
457            Traceback (most recent call last):
458            ...
459            ValueError: No 5th root of unity in self
460        """
461        if n == 1:
462            return self(1)
463        elif n == 2:
464            return self(-1)
465        raise ValueError, "No %sth root of unity in self"%n
466
467R = RealField()
468
469#*****************************************************************************
470#
471#     RealNumber -- element of Real Field
472#
473#     
474#
475#*****************************************************************************
476cdef class RealNumber(sage.structure.element.RingElement):
477    """
478    A real number.
479
480    Real numbers are printed to slightly less digits than their
481    internal precision, in order to avoid confusing roundoff issues
482    that occur because numbers are stored internally in binary.
483    """
484    cdef RealNumber _new(self):
485        """
486        Return a new real number with same parent as self.
487        """
488        cdef RealNumber x
489        x = PY_NEW(RealNumber)
490        x._parent = self._parent
491        mpfr_init2(x.value, (<RealField>self._parent).__prec)
492        x.init = 1
493        return x
494       
495    def __init__(self, RealField parent, x=0, int base=10):
496        """
497        Create a real number.  Should be called by first creating
498        a RealField, as illustrated in the examples.
499
500        EXAMPLES:
501            sage: R = RealField()
502            sage: R('1.2456')
503            1.24560000000000
504            sage: R = RealField(3)
505            sage: R('1.2456')
506            1.2
507
508        EXAMPLE: Rounding Modes
509            sage: w = RealField(3)(5/2)
510            sage: RealField(2, rnd="RNDZ")(w).str(2)
511            '10'
512            sage: RealField(2, rnd="RNDD")(w).str(2)
513            '10'
514            sage: RealField(2, rnd="RNDU")(w).str(2)
515            '11'
516            sage: RealField(2, rnd="RNDN")(w).str(2)
517            '10'
518
519        NOTES: A real number is an arbitrary precision mantissa with a
520        limited precision exponent.  A real number can have three
521        special values: Not-a-Number (NaN) or plus or minus
522        Infinity. NaN represents an uninitialized object, the result
523        of an invalid operation (like 0 divided by 0), or a value that
524        cannot be determined (like +Infinity minus
525        +Infinity). Moreover, like in the IEEE 754-1985 standard, zero
526        is signed, i.e. there are both +0 and -0; the behavior is the
527        same as in the IEEE 754-1985 standard and it is generalized to
528        the other functions supported by MPFR.
529       
530        """
531        self.init = 0
532        if parent is None:
533            raise TypeError
534        self._parent = parent
535        mpfr_init2(self.value, parent.__prec)
536        self.init = 1
537        if x is None: return
538        self._set(x, base)
539
540    cdef _set(self, x, int base):
541        # This should not be called except when the number is being created.
542        # Real Numbers are supposed to be immutable.
543        cdef RealNumber _x, n, d
544        cdef Integer _ix
545        cdef RealField parent
546        cdef gen _gen
547        parent = self._parent
548        if PY_TYPE_CHECK(x, RealNumber):
549            _x = x  # so we can get at x.value
550            mpfr_set(self.value, _x.value, parent.rnd)
551        elif PY_TYPE_CHECK(x, Integer):
552            mpfr_set_z(self.value, (<Integer>x).value, parent.rnd)           
553        elif PY_TYPE_CHECK(x, Rational):
554            mpfr_set_q(self.value, (<Rational>x).value, parent.rnd)
555        elif PY_TYPE_CHECK(x, gen) and x.type() == "t_REAL":
556            _gen = x
557            self._set_from_GEN_REAL(_gen.g)
558        elif isinstance(x, (int, long)):
559            _ix = Integer(x)
560            mpfr_set_z(self.value, _ix.value, parent.rnd)
561        #elif hasattr(x, '_mpfr_'):
562        #    return x._mpfr_(self)
563        else:
564            s = str(x).replace(' ','')
565            if mpfr_set_str(self.value, s, base, parent.rnd):
566                if s == 'NaN' or s == '@NaN@':
567                    mpfr_set_nan(self.value)
568                elif s == '+infinity':
569                    mpfr_set_inf(self.value, 1)
570                elif s == '-infinity':
571                    mpfr_set_inf(self.value, -1)
572                else:
573                    raise TypeError, "Unable to convert x (='%s') to real number."%s
574
575    cdef _set_from_GEN_REAL(self, GEN g):
576        """
577        EXAMPLES:
578            sage: rt2 = sqrt(pari('2.0'))
579            sage: rt2
580            1.414213562373095048801688724              # 32-bit
581            1.4142135623730950488016887242096980786    # 64-bit
582            sage: rt2.python()
583            1.414213562373095048801688724              # 32-bit
584            1.4142135623730950488016887242096980785    # 64-bit
585            sage: rt2.python().prec()
586            96                                         # 32-bit
587            128                                        # 64-bit
588            sage: pari(rt2.python()) == rt2
589            True
590            sage: for i in xrange(1, 1000):
591            ...       assert(sqrt(pari(i)) == pari(sqrt(pari(i)).python()))
592            sage: (-3.1415)._pari_().python()
593            -3.14150000000000018
594        """
595        cdef int sgn
596        sgn = signe(g)
597
598        if sgn == 0:
599            mpfr_set_ui(self.value, 0, GMP_RNDN)
600            return
601
602        cdef int wordsize
603
604        if sage.misc.misc.is_64_bit:
605            wordsize = 64
606        else:
607            wordsize = 32
608
609        cdef mpz_t mantissa
610        mpz_init(mantissa)
611
612        mpz_import(mantissa, lg(g) - 2, 1, wordsize/8, 0, 0, &g[2])
613
614        cdef int exponent
615        exponent = expo(g)
616
617        # Round to nearest for best results when setting a low-precision
618        # MPFR from a high-precision GEN
619        mpfr_set_z(self.value, mantissa, GMP_RNDN)
620        mpfr_mul_2si(self.value, self.value, exponent - wordsize * (lg(g) - 2) + 1, GMP_RNDN)
621
622        if sgn < 0:
623            mpfr_neg(self.value, self.value, GMP_RNDN)
624
625        mpz_clear(mantissa)
626
627    def __reduce__(self):
628        """
629        EXAMPLES:
630            sage: R = RealField(sci_not=1, prec=200, rnd='RNDU')
631            sage: b = R('393.39203845902384098234098230948209384028340')
632            sage: loads(dumps(b)) == b
633            True
634            sage: b = R(1)/R(0); b
635            +infinity
636            sage: loads(dumps(b)) == b
637            True
638            sage: b = R(-1)/R(0); b
639            -infinity
640            sage: loads(dumps(b)) == b
641            True
642            sage: b = R(-1).sqrt(); b
643            1.0000000000000000000000000000000000000000000000000000000000*I
644            sage: loads(dumps(b)) == b
645            True
646        """
647        s = self.str(32, no_sci=False, e='@')
648        return (__create__RealNumber_version0, (self._parent, s, 32))
649
650    def  __dealloc__(self):
651        if self.init:
652            mpfr_clear(self.value)
653
654    def __repr__(self):
655        return self.str(10)
656
657    def _latex_(self):
658        return str(self)
659
660    def _interface_init_(self):
661        """
662        Return string representation of self in base 10, avoiding
663        scientific notation except for very large or very small numbers.
664
665        This is most likely to make sense in other computer algebra
666        systems (this function is the default for exporting to other
667        computer algebra systems).
668
669        EXAMPLES:
670            sage: n = 1.3939494594
671            sage: n._interface_init_()
672            '1.39394945939999'
673        """
674        return self.str(10, no_sci=True)
675
676    def __hash__(self):
677        return hash(self.str(16))
678
679    def _im_gens_(self, codomain, im_gens):
680        return codomain(self) # since 1 |--> 1
681
682    def real(self):
683        """
684        Return the real part of self.
685
686        (Since self is a real number, this simply returns self.)
687        """
688        return self
689
690    def parent(self):
691        """
692        EXAMPLES:
693            sage: R = RealField()
694            sage: a = R('1.2456')
695            sage: a.parent()
696            Real Field with 53 bits of precision
697        """
698        return self._parent
699
700    def str(self, int base=10, no_sci=None, e=None, int truncate=1):
701        """
702        INPUT:
703             base -- base for output
704             no_sci -- if 2, never print using scientific notation;
705                       if 1 or True, print using scientific notation only
706                       for very large or very small numbers;
707                       if 0 or False always print with scientific notation;
708                       if None (the default), print how the parent prints.
709             e - symbol used in scientific notation; defaults to 'e' for
710                       base<=10, and '@' otherwise
711             truncate -- if True, round off the last digits in printing to
712                       lessen confusing base-2 roundoff issues.
713                         
714        EXAMPLES:
715            sage: a = 61/3.0; a
716            20.3333333333333
717            sage: a.str(truncate=False)
718            '20.333333333333332'
719            sage: a.str(2)
720            '10100.010101010101010101010101010101010101010101010101'
721            sage: a.str(no_sci=False)
722            '2.03333333333333e1'       
723            sage: a.str(16, no_sci=False)
724            '1.4555555555555@1'
725            sage: b = 2.0^99
726            sage: b.str()
727            '633825300114114000000000000000'
728            sage: b.str(no_sci=False)
729            '6.33825300114114e29'
730            sage: b.str(no_sci=True)
731            '633825300114114000000000000000'
732            sage: c = 2.0^100
733            sage: c.str()
734            '1.26765060022822e30'
735            sage: c.str(no_sci=False)
736            '1.26765060022822e30'
737            sage: c.str(no_sci=True)
738            '1.26765060022822e30'
739            sage: c.str(no_sci=2)
740            '1267650600228220000000000000000'
741            sage: 0.5^53
742            0.000000000000000111022302462515
743            sage: 0.5^54
744            5.55111512312578e-17
745        """
746        if base < 2 or base > 36:
747            raise ValueError, "the base (=%s) must be between 2 and 36"%base
748        if mpfr_nan_p(self.value):
749            if base >= 24:
750                return "@NaN@"
751            else:
752                return "NaN"
753        elif mpfr_inf_p(self.value):
754            if mpfr_sgn(self.value) > 0:
755                return "+infinity"
756            else:
757                return "-infinity"
758       
759        if e is None:
760            if base > 10:
761                e = '@'
762            else:
763                e = 'e'
764
765        cdef char *s
766        cdef mp_exp_t exponent
767
768        cdef int reqdigits
769
770        reqdigits = 0
771
772        if base == 10 and truncate:
773
774            # This computes reqdigits == floor(log_{10}(2^(b-1))),
775            # which is the number of *decimal* digits that are
776            # "right", given that the last binary bit of the binary
777            # number can be off.  That is, if this real is within a
778            # relative error of 2^(-b) of an exact decimal with
779            # reqdigits digits, that decimal will be returned.
780            # This is equivalent to saying that exact decimals with
781            # reqdigits digits differ by at least 2*2^(-b) (relative).
782
783            # (Depending on the precision and the exact number involved,
784            # adjacent exact decimals can differ by far more than 2*2^(-b)
785            # (relative).)
786
787            # This avoids the confusion a lot of people have with the last
788            # 1-2 binary digits being wrong due to rounding coming from
789            # representating numbers in binary.
790
791            reqdigits = ((<RealField>self._parent).__prec - 1) * 0.3010299956
792            if reqdigits <= 1: reqdigits = 2
793       
794        _sig_on
795        s = mpfr_get_str(<char*>0, &exponent, base, reqdigits,
796                         self.value, (<RealField>self._parent).rnd)
797        _sig_off
798        if s == <char*> 0:
799            raise RuntimeError, "Unable to convert an mpfr number to a string."
800        t = str(s)
801        free(s)
802       
803       
804        cdef int digits
805        digits = len(t)
806
807        if no_sci is None:
808            no_sci = not (<RealField>self._parent).sci_not
809
810        if no_sci==True and (-exponent > digits or exponent > 2*digits):
811            no_sci = False
812
813        if no_sci==False:
814            if t[0] == "-":
815                return "-%s.%s%s%s"%(t[1:2], t[2:], e, exponent-1)
816            return "%s.%s%s%s"%(t[0], t[1:], e, exponent-1)
817
818        lpad = ''
819        if exponent <= 0:
820            n = len(t)
821            lpad = '0.' + '0'*abs(exponent)
822        else:
823            n = exponent
824        if t[0] == '-':
825            lpad = '-' + lpad
826            t = t[1:]
827        z = lpad + str(t[:n])
828        w = t[n:]
829        if len(w) > 0:
830            z = z + ".%s"%w
831        elif exponent > 0:
832            z = z + '0'*(n-len(t))
833        return z
834
835    def __copy__(self):
836        """
837        Return copy of self -- since self is immutable, we just return self again.
838
839        EXAMPLES:
840            sage: a = 3.5
841            sage: copy(a) is  a
842            True
843        """
844        return self    # since object is immutable.
845
846    def integer_part(self):
847        """
848        If in decimal this number is written n.defg, returns n.
849
850        OUTPUT:
851            -- a SAGE Integer
852
853        EXAMPLE:
854            sage: a = 119.41212
855            sage: a.integer_part()
856            119
857
858        A big number with no decimal point:
859            sage: a = RR(10^17); a
860            100000000000000000
861            sage: a.integer_part()
862            100000000000000000
863        """
864        s = self.str(base=32, no_sci=True)
865        i = s.find(".")
866        if i != -1:
867            return Integer(s[:i], base=32)
868        else:
869            return Integer(s, base=32)
870
871    ########################
872    #   Basic Arithmetic
873    ########################
874
875    cdef ModuleElement _add_c_impl(self, ModuleElement other):
876        """
877        Add two real numbers with the same parent.
878       
879        EXAMPLES:
880            sage: R = RealField()
881            sage: R(-1.5) + R(2.5)
882            1.00000000000000
883        """
884        cdef RealNumber x
885        x = self._new()
886        mpfr_add(x.value, self.value, (<RealNumber>other).value, (<RealField>self._parent).rnd)
887        return x
888       
889    def __invert__(self):
890        return self._parent(1) / self
891   
892    cdef ModuleElement _sub_c_impl(self, ModuleElement right):
893        """
894        Subtract two real numbers with the same parent.
895       
896        EXAMPLES:
897            sage: R = RealField()
898            sage: R(-1.5) - R(2.5)
899            -4.00000000000000
900        """
901        cdef RealNumber x
902        x = self._new()
903        mpfr_sub(x.value, self.value, (<RealNumber>right).value, (<RealField> self._parent).rnd)
904        return x
905       
906    cdef RingElement _mul_c_impl(self, RingElement right):
907        """
908        Multiply two real numbers with the same parent.
909       
910        EXAMPLES:
911            sage: R = RealField()
912            sage: R(-1.5) * R(2.5)
913            -3.75000000000000
914
915        If two elements have different precision, arithmetic
916        operations are performed after coercing to the lower
917        precision.
918       
919            sage: R20 = RealField(20)
920            sage: R100 = RealField(100)
921            sage: a = R20('393.3902834028345')
922            sage: b = R100('393.3902834028345')
923            sage: a
924            393.39
925            sage: b
926            393.39028340283450000000000000
927            sage: a*b
928            154750
929            sage: b*a
930            154750
931            sage: parent(b*a)
932            Real Field with 20 bits of precision
933        """
934        cdef RealNumber x
935        x = self._new()
936        mpfr_mul(x.value, self.value, (<RealNumber>right).value, (<RealField>self._parent).rnd)
937        return x
938   
939       
940    cdef RingElement _div_c_impl(self, RingElement right):
941        """
942        Divide self by other, where both are real numbers with the same parent.
943
944        EXAMPLES:
945            sage: RR(1)/RR(3)
946            0.333333333333333
947            sage: RR(1)/RR(0)
948            +infinity       
949
950            sage: R = RealField()
951            sage: R(-1.5) / R(2.5)
952            -0.599999999999999
953        """
954        cdef RealNumber x
955        x = self._new()
956        mpfr_div((<RealNumber>x).value, self.value,
957                 (<RealNumber>right).value, (<RealField>self._parent).rnd)
958        return x
959   
960    cdef ModuleElement _neg_c_impl(self):
961        cdef RealNumber x
962        x = self._new()
963        mpfr_neg(x.value, self.value, (<RealField>self._parent).rnd)
964        return x
965
966    def __abs__(self):
967        return self.abs()
968
969    cdef RealNumber abs(RealNumber self):
970        cdef RealNumber x
971        x = self._new()       
972        mpfr_abs(x.value, self.value, (<RealField>self._parent).rnd)
973        return x
974
975    # Bit shifting
976    def _lshift_(RealNumber self, n):
977        cdef RealNumber x
978        if n > sys.maxint:
979            raise OverflowError, "n (=%s) must be <= %s"%(n, sys.maxint)
980        x = self._new()
981        mpfr_mul_2exp(x.value, self.value, n, (<RealField>self._parent).rnd)
982        return x
983
984    def __lshift__(x, y):
985        """
986        EXAMPLES:
987            sage: 1.0 << 32
988            4294967296.00000
989        """
990        if isinstance(x, RealNumber) and isinstance(y, (int,long, Integer)):
991            return x._lshift_(y)
992        return sage.structure.coerce.bin_op(x, y, operator.lshift)
993
994    def _rshift_(RealNumber self, n):
995        if n > sys.maxint:
996            raise OverflowError, "n (=%s) must be <= %s"%(n, sys.maxint)
997        cdef RealNumber x
998        x = self._new()
999        mpfr_div_2exp(x.value, self.value, n, (<RealField>self._parent).rnd)
1000        return x
1001
1002    def __rshift__(x, y):
1003        """
1004        EXAMPLES:
1005            sage: 1024.0 >> 7
1006            8.00000000000000
1007        """
1008        if isinstance(x, RealNumber) and isinstance(y, (int,long,Integer)):
1009            return x._rshift_(y)
1010        return sage.structure.coerce.bin_op(x, y, operator.rshift)
1011   
1012    def multiplicative_order(self):
1013        if self == 1:
1014            return 1
1015        elif self == -1:
1016            return -1
1017        return sage.rings.infinity.infinity
1018
1019    def sign(self):
1020        return mpfr_sgn(self.value)
1021
1022    def prec(self):
1023        return (<RealField>self._parent).__prec
1024
1025    ###################
1026    # Rounding etc
1027    ###################
1028
1029    def round(self):
1030        """
1031        Rounds self to the nearest real number. There are 4
1032        rounding modes. They are
1033
1034        EXAMPLES:
1035            RNDN -- round to nearest:
1036           
1037            sage: R = RealField(20,False,'RNDN')
1038            sage: R(22.454)
1039            22.454
1040            sage: R(22.491)
1041            22.490
1042
1043            RNDZ -- round towards zero:
1044            sage: R = RealField(20,False,'RNDZ')
1045            sage: R(22.454)
1046            22.453
1047            sage: R(22.491)
1048            22.490
1049                       
1050            RNDU -- round towards plus infinity:
1051            sage: R = RealField(20,False,'RNDU')
1052            sage: R(22.454)
1053            22.454
1054            sage: R(22.491)
1055            22.491
1056
1057            RNDU -- round towards minus infinity:
1058            sage: R = RealField(20,False,'RNDD')
1059            sage: R(22.454)
1060            22.453
1061            sage: R(22.491)
1062            22.490
1063        """
1064        cdef RealNumber x
1065        x = self._new()
1066        mpfr_round(x.value, self.value)
1067        return x
1068
1069    def floor(self):
1070        """
1071        Returns the floor of this number
1072
1073        EXAMPLES:
1074            sage: R = RealField()
1075            sage: (2.99).floor()
1076            2
1077            sage: (2.00).floor()
1078            2
1079            sage: floor(RR(-5/2))
1080            -3
1081        """
1082        cdef RealNumber x
1083        x = self._new()
1084        mpfr_floor(x.value, self.value)
1085        return x.integer_part()
1086
1087    def ceil(self):
1088        """
1089        Returns the ceiling of this number
1090
1091        OUTPUT:
1092            integer
1093
1094        EXAMPLES:
1095            sage: (2.99).ceil()
1096            3
1097            sage: (2.00).ceil()
1098            2
1099            sage: (2.01).ceil()
1100            3
1101
1102            sage: ceil(10^16 * 1.0)
1103            10000000000000000
1104            sage: ceil(10^17 * 1.0)
1105            100000000000000000
1106        """
1107        cdef RealNumber x
1108        x = self._new()
1109        mpfr_ceil(x.value, self.value)
1110        return x.integer_part()
1111
1112    def ceiling(self):
1113        return self.ceil()
1114
1115    def trunc(self):
1116        """
1117        Truncates this number
1118
1119        EXAMPLES:
1120            sage: (2.99).trunc()
1121            2.00000000000000
1122            sage: (-0.00).trunc()
1123            -0.000000000000000
1124            sage: (0.00).trunc()
1125            0.000000000000000
1126        """
1127        cdef RealNumber x
1128        x = self._new()
1129        mpfr_trunc(x.value, self.value)
1130        return x
1131
1132    def frac(self):
1133        """
1134        frac returns a real number > -1 and < 1. it satisfies the
1135        relation:
1136            x = x.trunc() + x.frac()
1137
1138        EXAMPLES:
1139            sage: (2.99).frac()
1140            0.990000000000000
1141            sage: (2.50).frac()
1142            0.500000000000000
1143            sage: (-2.79).frac()
1144            -0.790000000000000
1145        """
1146        cdef RealNumber x
1147        x = self._new()
1148        mpfr_frac(x.value, self.value, (<RealField>self._parent).rnd)
1149        return x
1150
1151    ###########################################
1152    # Conversions
1153    ###########################################
1154
1155    def __float__(self):
1156        return mpfr_get_d(self.value, (<RealField>self._parent).rnd)
1157   
1158    def __int__(self):
1159        """
1160        Returns integer truncation of this real number.
1161        """
1162        s = self.str(32)
1163        i = s.find('.')
1164        return int(s[:i], 32)
1165
1166    def __long__(self):
1167        """
1168        Returns long integer truncation of this real number.
1169        """
1170        s = self.str(32)
1171        i = s.find('.')
1172        return long(s[:i], 32)
1173
1174    def __complex__(self):
1175        return complex(float(self))
1176
1177    def _complex_number_(self):
1178        return sage.rings.complex_field.ComplexField(self.prec())(self)
1179
1180    def _pari_(self):
1181        """
1182        Returns self as a Pari floating-point number.
1183
1184        EXAMPLES:
1185            sage: RR(2.0)._pari_()
1186            2.000000000000000000
1187            sage: RealField(250).pi()._pari_()
1188            3.141592653589793238462643383
1189            sage: RR(0.0)._pari_()
1190            0.E-19
1191            sage: RR(-1.234567)._pari_()
1192            -1.2345669999999999700
1193            sage: RR(2.0).sqrt()._pari_()
1194            1.4142135623730951455
1195            sage: RR(2.0).sqrt()._pari_().python()
1196            1.41421356237309514
1197            sage: RR(2.0).sqrt()._pari_().python().prec()
1198            64
1199            sage: RealField(70)(pi)._pari_().python().prec()
1200            96                                         # 32-bit
1201            128                                        # 64-bit
1202            sage: for i in xrange(1, 1000):
1203            ...       assert(RR(i).sqrt() == RR(i).sqrt()._pari_().python())
1204        """
1205        # return sage.libs.pari.all.pari.new_with_bits_prec(str(self), (<RealField>self._parent).__prec)
1206
1207        # This uses interfaces of MPFR and Pari which are documented
1208        # (and not marked subject-to-change).  It could be faster
1209        # by using internal interfaces of MPFR, which are documented
1210        # as subject-to-change.
1211
1212        if mpfr_nan_p(self.value) or mpfr_inf_p(self.value):
1213            raise ValueError, 'Cannot convert NaN or infinity to Pari float'
1214
1215        cdef int wordsize
1216
1217        if sage.misc.misc.is_64_bit:
1218            wordsize = 64
1219        else:
1220            wordsize = 32
1221
1222        cdef int prec
1223        prec = (<RealField>self._parent).__prec
1224
1225        # We round up the precision to the nearest multiple of wordsize.
1226        cdef int rounded_prec
1227        rounded_prec = (self.prec() + wordsize - 1) & ~(wordsize - 1)
1228
1229        # Yes, assigning to self works fine, even in Pyrex.
1230        if rounded_prec > prec:
1231            self = RealField(rounded_prec)(self)
1232
1233        # Now we can extract the mantissa, and it will be normalized
1234        # (the most significant bit of the most significant word will be 1).
1235        cdef mpz_t mantissa
1236        cdef mp_exp_t exponent
1237        mpz_init(mantissa)
1238
1239        exponent = mpfr_get_z_exp(mantissa, self.value)
1240
1241        cdef GEN pari_float
1242        pari_float = cgetr(2 + rounded_prec / wordsize)
1243
1244        mpz_export(&pari_float[2], NULL, 1, wordsize/8, 0, 0, mantissa)
1245
1246        if mpfr_zero_p(self.value):
1247            setexpo(pari_float, -rounded_prec)
1248        else:
1249            setexpo(pari_float, exponent + rounded_prec - 1)
1250        setsigne(pari_float, mpfr_sgn(self.value))
1251
1252        cdef PariInstance P
1253        P = sage.libs.pari.all.pari
1254       
1255        gen = P.new_gen(pari_float)
1256
1257        mpz_clear(mantissa)
1258
1259        return gen
1260
1261    ###########################################
1262    # Comparisons: ==, !=, <, <=, >, >=
1263    ###########################################
1264
1265    def is_NaN(self):
1266        return bool(mpfr_nan_p(self.value))
1267
1268    def is_positive_infinity(self):
1269        """
1270        EXAMPLES:
1271            sage: a = RR('1.494') / RR(0); a
1272            +infinity
1273            sage: a.is_positive_infinity()
1274            True
1275            sage: a = -RR('1.494') / RR(0); a
1276            -infinity
1277            sage: RR(1.5).is_positive_infinity()
1278            False
1279            sage: a.is_positive_infinity()
1280            False
1281        """
1282        return bool(mpfr_inf_p(self.value) and mpfr_sgn(self.value) > 0)
1283
1284    def is_negative_infinity(self):
1285        """
1286        EXAMPLES:
1287            sage: a = RR('1.494') / RR(0); a
1288            +infinity
1289            sage: a.is_negative_infinity()
1290            False
1291            sage: a = -RR('1.494') / RR(0); a
1292            -infinity
1293            sage: RR(1.5).is_negative_infinity()
1294            False
1295            sage: a.is_negative_infinity()
1296            True
1297        """
1298        return bool(mpfr_inf_p(self.value) and mpfr_sgn(self.value) < 0)       
1299
1300    def is_infinity(self):
1301        """
1302        EXAMPLES:
1303            sage: a = RR('1.494') / RR(0); a
1304            +infinity
1305            sage: a.is_infinity()
1306            True
1307            sage: a = -RR('1.494') / RR(0); a
1308            -infinity
1309            sage: a.is_infinity()
1310            True
1311            sage: RR(1.5).is_infinity()
1312            False
1313        """
1314        return bool(mpfr_inf_p(self.value))
1315
1316    def __richcmp__(left, right, int op):
1317        return (<RingElement>left)._richcmp(right, op)
1318
1319    cdef int _cmp_c_impl(left, Element right) except -2:
1320        cdef RealNumber self, x
1321        self = left
1322        x = right
1323
1324        cdef int a,b
1325        a = mpfr_nan_p(self.value)
1326        b = mpfr_nan_p(x.value)
1327        if a != b:
1328            return -1    # nothing is equal to Nan
1329        cdef int i
1330        i = mpfr_cmp(self.value, x.value)
1331        if i < 0:
1332            return -1
1333        elif i == 0:
1334            return 0
1335        else:
1336            return 1
1337   
1338
1339    ############################
1340    # Special Functions
1341    ############################
1342
1343    def sqrt(self):
1344        """
1345        Return a square root of self.
1346
1347        If self is negative a complex number is returned.
1348
1349        If you use self.square_root() then a real number will always
1350        be returned (though it will be NaN if self is negative).
1351
1352        EXAMPLES:
1353            sage: r = 4.0
1354            sage: r.sqrt()
1355            2.00000000000000
1356            sage: r.sqrt()^2 == r
1357            True
1358
1359            sage: r = 4344
1360            sage: r.sqrt()
1361            65.9090282131363
1362            sage: r.sqrt()^2 == r
1363            True
1364            sage: r.sqrt()^2 - r           
1365             0.000000000000000
1366
1367            sage: r = -2.0
1368            sage: r.sqrt()
1369            1.41421356237309*I
1370            """
1371        if self >= 0:
1372            return self.square_root()
1373        return self._complex_number_().sqrt()
1374       
1375
1376    def square_root(self):
1377        """
1378        Return a square root of self.  A real number will always be
1379        returned (though it will be NaN if self is negative).
1380
1381        Use self.sqrt() to get a complex number if self is negative.
1382
1383        EXAMPLES:
1384            sage: r = -2.0
1385            sage: r.square_root()
1386            NaN
1387            sage: r.sqrt()
1388            1.41421356237309*I
1389        """
1390        cdef RealNumber x
1391        x = self._new()
1392        _sig_on               
1393        mpfr_sqrt(x.value, self.value, (<RealField>self._parent).rnd)
1394        _sig_off
1395        return x
1396
1397    def cube_root(self):
1398        """
1399        Return the cubic root (defined over the real numbers) of self.
1400
1401        EXAMPLES:
1402            sage: r = 125.0; r.cube_root()
1403            5.00000000000000
1404            sage: r = -119.0
1405            sage: r.cube_root()^3 - r       # illustrates precision loss
1406            -0.0000000000000142108547152020
1407        """
1408        cdef RealNumber x
1409        x = self._new()
1410        _sig_on               
1411        mpfr_cbrt(x.value, self.value, (<RealField>self._parent).rnd)
1412        _sig_off       
1413        return x
1414
1415    def __pow(self, RealNumber exponent):
1416        cdef RealNumber x
1417        x = self._new()
1418        _sig_on               
1419        mpfr_pow(x.value, self.value, exponent.value, (<RealField>self._parent).rnd)
1420        _sig_off
1421        if mpfr_nan_p(x.value):
1422            return self._complex_number_()**exponent._complex_number_()
1423        return x
1424
1425    def __pow__(self, exponent, modulus):
1426        """
1427        Compute self raised to the power of exponent, rounded in
1428        the direction specified by the parent of self.
1429
1430        If the result is not a real number, self and the exponent are
1431        both coerced to complex numbers (with sufficient precision),
1432        then the exponentiation is computed in the complex numbers.
1433        Thus this function can return either a real or complex number.
1434       
1435        EXAMPLES:
1436            sage: R = RealField(30)
1437            sage: a = R('1.23456')
1438            sage: a^20
1439            67.646297
1440            sage: a^a
1441            1.2971114
1442            sage: b = R(-1)
1443            sage: b^(1/2)     
1444            1.0000000*I                   # 32-bit
1445            -0.00000000000000000010842021 + 1.0000000*I   # 64-bit
1446        """
1447        cdef RealNumber x
1448        if not PY_TYPE_CHECK(self, RealNumber):
1449            return self.__pow__(float(exponent))
1450        if not PY_TYPE_CHECK(exponent, RealNumber):
1451            x = self
1452            exponent = x._parent(exponent)
1453        return self.__pow(exponent)
1454
1455    def log(self, base='e'):
1456        """
1457        EXAMPLES:
1458            sage: R = RealField()
1459            sage: R(2).log()
1460            0.693147180559945
1461        """
1462        cdef RealNumber x
1463        if base == 'e':
1464            x = self._new()
1465            _sig_on               
1466            mpfr_log(x.value, self.value, (<RealField>self._parent).rnd)
1467            _sig_off       
1468            return x
1469        elif base == 10:
1470            return self.log10()
1471        elif base == 2:
1472            return self.log2()
1473        else:
1474            return self.log() / (self.parent()(base)).log()
1475
1476    def log2(self):
1477        """
1478        Returns log to the base 2 of self
1479
1480        EXAMPLES:
1481            sage: r = 16.0
1482            sage: r.log2()
1483            4.00000000000000
1484
1485            sage: r = 31.9; r.log2()
1486            4.99548451887750
1487
1488            sage: r = 0.0
1489            sage: r.log2()
1490            -infinity
1491        """
1492        cdef RealNumber x
1493        x = self._new()
1494        _sig_on               
1495        mpfr_log2(x.value, self.value, (<RealField>self._parent).rnd)
1496        _sig_off       
1497        return x
1498
1499    def log10(self):
1500        """
1501        Returns log to the base 10 of self
1502
1503        EXAMPLES:
1504            sage: r = 16.0; r.log10()
1505            1.20411998265592
1506            sage: r.log() / log(10)
1507            1.20411998265804
1508
1509            sage: r = 39.9; r.log10()
1510            1.60097289568674
1511
1512            sage: r = 0.0
1513            sage: r.log10()
1514            -infinity
1515
1516            sage: r = -1.0
1517            sage: r.log10()
1518            NaN
1519
1520        """
1521        cdef RealNumber x
1522        x = self._new()
1523        _sig_on               
1524        mpfr_log10(x.value, self.value, (<RealField>self._parent).rnd)
1525        _sig_off
1526        return x
1527
1528    def exp(self):
1529        r"""
1530        Returns $e^\code{self}$
1531
1532        EXAMPLES:
1533            sage: r = 0.0
1534            sage: r.exp()
1535            1.00000000000000
1536
1537            sage: r = 32.3
1538            sage: a = r.exp(); a
1539            106588847274864
1540            sage: a.log()
1541            32.2999999999999
1542
1543            sage: r = -32.3
1544            sage: r.exp()
1545            0.00000000000000938184458849868
1546        """
1547        cdef RealNumber x
1548        x = self._new()
1549        _sig_on               
1550        mpfr_exp(x.value, self.value, (<RealField>self._parent).rnd)
1551        _sig_off       
1552        return x
1553
1554    def exp2(self):
1555        """
1556        Returns $2^\code{self}$
1557
1558        EXAMPLES:
1559            sage: r = 0.0
1560            sage: r.exp2()
1561            1.00000000000000
1562           
1563            sage: r = 32.0
1564            sage: r.exp2()
1565            4294967296.00000
1566
1567            sage: r = -32.3
1568            sage: r.exp2()
1569            0.000000000189117248253020
1570
1571        """
1572        cdef RealNumber x
1573        x = self._new()
1574        _sig_on               
1575        mpfr_exp2(x.value, self.value, (<RealField>self._parent).rnd)
1576        _sig_off       
1577        return x
1578
1579    def exp10(self):
1580        r"""
1581        Returns $10^\code{self}$
1582
1583        EXAMPLES:
1584            sage: r = 0.0
1585            sage: r.exp10()
1586            1.00000000000000
1587           
1588            sage: r = 32.0
1589            sage: r.exp10()
1590            1.00000000000000e32
1591
1592            sage: r = -32.3
1593            sage: r.exp10()
1594            5.01187233627275e-33
1595        """
1596        cdef RealNumber x
1597        x = self._new()
1598        _sig_on               
1599        mpfr_exp10(x.value, self.value, (<RealField>self._parent).rnd)
1600        _sig_off       
1601        return x
1602
1603    def cos(self):
1604        """
1605        Returns the cosine of this number
1606
1607        EXAMPLES:
1608            sage: t=RR.pi()/2
1609            sage: t.cos()
1610            6.12323399573676e-17
1611        """
1612        cdef RealNumber x
1613        x = self._new()
1614        _sig_on               
1615        mpfr_cos(x.value, self.value, (<RealField>self._parent).rnd)
1616        _sig_off       
1617        return x
1618
1619    ##########################################################
1620    # it would be nice to get zero back here:
1621    # sage: R(-1).acos().sin()
1622    # _57 = -0.50165576126683320234e-19
1623    # i think this could be "fixed" by using MPFI. (put on to-do list.)
1624    #
1625    # this seems to work ok:
1626    # sage: R(-1).acos().cos()
1627    # _58 = -0.10000000000000000000e1
1628    def sin(self):
1629        """
1630        Returns the sine of this number
1631       
1632        EXAMPLES:
1633            sage: R = RealField(100)
1634            sage: R(2).sin()
1635            0.90929742682568169539601986591
1636        """
1637        cdef RealNumber x
1638        x = self._new()
1639        _sig_on               
1640        mpfr_sin(x.value, self.value, (<RealField>self._parent).rnd)
1641        _sig_off
1642        return x
1643   
1644    def tan(self):
1645        """
1646        Returns the tangent of this number
1647       
1648        EXAMPLES:
1649            sage: q = RR.pi()/3
1650            sage: q.tan()
1651            1.73205080756887
1652            sage: q = RR.pi()/6
1653            sage: q.tan()
1654            0.577350269189625
1655        """
1656        cdef RealNumber x
1657        x = self._new()
1658        _sig_on               
1659        mpfr_tan(x.value, self.value, (<RealField>self._parent).rnd)
1660        _sig_off
1661        return x
1662
1663    def sincos(self):
1664        """
1665        Returns a pair consisting of the sine and cosine.
1666
1667        EXAMPLES:
1668            sage: R = RealField()
1669            sage: t = R.pi()/6
1670            sage: t.sincos()
1671            (0.499999999999999, 0.866025403784438)
1672        """
1673        cdef RealNumber x,y
1674        x = self._new()
1675        y = self._new()
1676        _sig_on               
1677        mpfr_sin_cos(x.value, y.value, self.value, (<RealField>self._parent).rnd)
1678        _sig_off       
1679        return x,y
1680
1681
1682    # int mpfr_sin_cos (mpfr_t rop, mpfr_t op, mpfr_t, mp_rnd_t rnd)
1683
1684    def acos(self):
1685        """
1686        Returns the inverse cosine of this number
1687
1688        EXAMPLES:
1689            sage: q = RR.pi()/3
1690            sage: i = q.cos()
1691            sage: i.acos() == q
1692            True
1693        """
1694        cdef RealNumber x
1695        x = self._new()
1696        _sig_on               
1697        mpfr_acos(x.value, self.value, (<RealField>self._parent).rnd)
1698        _sig_off       
1699        return x
1700
1701    def asin(self):
1702        """
1703        Returns the inverse sine of this number
1704
1705        EXAMPLES:
1706            sage: q = RR.pi()/5
1707            sage: i = q.sin()
1708            sage: i.asin() == q
1709            True
1710            sage: i.asin() - q
1711            0.000000000000000
1712        """
1713        cdef RealNumber x
1714        x = self._new()
1715        _sig_on               
1716        mpfr_asin(x.value, self.value, (<RealField>self._parent).rnd)
1717        _sig_off       
1718        return x
1719
1720    def atan(self):
1721        """
1722        Returns the inverse tangent of this number
1723
1724        EXAMPLES:
1725            sage: q = RR.pi()/5
1726            sage: i = q.tan()
1727        """
1728        cdef RealNumber x
1729        x = self._new()
1730        _sig_on               
1731        mpfr_atan(x.value, self.value, (<RealField>self._parent).rnd)
1732        _sig_off       
1733        return x
1734
1735    #int mpfr_acos _PROTO ((mpfr_ptr, mpfr_srcptr, mp_rnd_t));
1736    #int mpfr_asin _PROTO ((mpfr_ptr, mpfr_srcptr, mp_rnd_t));
1737    #int mpfr_atan _PROTO ((mpfr_ptr, mpfr_srcptr, mp_rnd_t));
1738
1739    def cosh(self):
1740        """
1741        Returns the hyperbolic cosine of this number
1742
1743        EXAMPLES:
1744            sage: q = RR.pi()/12
1745            sage: q.cosh()
1746            1.03446564009551
1747        """
1748        cdef RealNumber x
1749        x = self._new()
1750        _sig_on               
1751        mpfr_cosh(x.value, self.value, (<RealField>self._parent).rnd)
1752        _sig_off       
1753        return x
1754
1755    def sinh(self):
1756        """
1757        Returns the hyperbolic sine of this number
1758
1759        EXAMPLES:
1760            sage: q = RR.pi()/12
1761            sage: q.sinh()
1762            0.264800227602270
1763        """
1764        cdef RealNumber x
1765        x = self._new()
1766        _sig_on               
1767        mpfr_sinh(x.value, self.value, (<RealField>self._parent).rnd)
1768        _sig_off       
1769        return x
1770
1771    def tanh(self):
1772        """
1773        Returns the hyperbolic tangent of this number
1774
1775        EXAMPLES:
1776            sage: q = RR.pi()/11
1777            sage: q.tanh()
1778            0.278079429295850
1779        """
1780        cdef RealNumber x
1781        x = self._new()
1782        _sig_on               
1783        mpfr_tanh(x.value, self.value, (<RealField>self._parent).rnd)
1784        _sig_off       
1785        return x
1786
1787    def acosh(self):
1788        """
1789        Returns the hyperbolic inverse cosine of this number
1790
1791        EXAMPLES:
1792            sage: q = RR.pi()/2
1793            sage: i = q.cosh() ; i
1794            2.50917847865805
1795        """
1796        cdef RealNumber x
1797        x = self._new()
1798        _sig_on               
1799        mpfr_acosh(x.value, self.value, (<RealField>self._parent).rnd)
1800        _sig_off       
1801        return x
1802
1803    def asinh(self):
1804        """
1805        Returns the hyperbolic inverse sine of this number
1806
1807        EXAMPLES:
1808            sage: q = RR.pi()/7
1809            sage: i = q.sinh() ; i
1810            0.464017630492990
1811            sage: i.asinh() - q
1812            0.000000000000000
1813        """
1814        cdef RealNumber x
1815        x = self._new()
1816        _sig_on               
1817        mpfr_asinh(x.value, self.value, (<RealField>self._parent).rnd)
1818        _sig_off       
1819        return x
1820
1821    def atanh(self):
1822        """
1823        Returns the hyperbolic inverse tangent of this number
1824
1825        EXAMPLES:
1826            sage: q = RR.pi()/7
1827            sage: i = q.tanh() ; i
1828            0.420911241048534
1829            sage: i.atanh() - q
1830            0.000000000000000
1831        """
1832        cdef RealNumber x
1833        x = self._new()
1834        _sig_on       
1835        mpfr_atanh(x.value, self.value, (<RealField>self._parent).rnd)
1836        _sig_off       
1837        return x
1838
1839    def agm(self, other):
1840        """
1841        Return the arithmetic-geometric mean of self and other. The
1842        arithmetic-geometric mean is the common limit of the sequences
1843        $u_n$ and $v_n$, where $u_0$ is self, $v_0$ is other,
1844        $u_{n+1}$ is the arithmetic mean of $u_n$ and $v_n$, and
1845        $v_{n+1}$ is the geometric mean of u_n and v_n. If any operand
1846        is negative, the return value is \code{NaN}.
1847        """
1848        cdef RealNumber x, _other
1849        if not isinstance(other, RealNumber) or other.parent() != self._parent:
1850            _other = self._parent(other)
1851        else:
1852            _other = other
1853        x = self._new()
1854        _sig_on       
1855        mpfr_agm(x.value, self.value, _other.value, (<RealField>self._parent).rnd)
1856        _sig_off
1857        return x
1858
1859   
1860    def erf(self):
1861        """
1862        Returns the value of the error function on self.
1863
1864        EXAMPLES:                                                           
1865           sage: R = RealField()
1866           sage: R(6).erf()
1867           1.00000000000000
1868        """
1869        cdef RealNumber x
1870        x = self._new()
1871        _sig_on
1872        mpfr_erf(x.value, self.value, (<RealField>self._parent).rnd)
1873        _sig_off
1874        return x
1875
1876
1877    def gamma(self):
1878        """
1879        The Euler gamma function. Return gamma of self.
1880
1881        EXAMPLES:                                                           
1882           sage: R = RealField()                                               
1883           sage: R(6).gamma()
1884           120.000000000000
1885           sage: R(1.5).gamma()
1886           0.886226925452758
1887        """
1888        cdef RealNumber x
1889        x = self._new()
1890        _sig_on
1891        mpfr_gamma(x.value, self.value, (<RealField>self._parent).rnd)
1892        _sig_off
1893        return x
1894
1895    def zeta(self):
1896        r"""
1897        Return the Riemann zeta function evaluated at this real number.
1898
1899        \note{PARI is vastly more efficient at computing the Riemann zeta
1900        function.   See the example below for how to use it.}
1901
1902        EXAMPLES:                                                           
1903            sage: R = RealField()
1904            sage: R(2).zeta()
1905            1.64493406684822
1906            sage: R.pi()^2/6
1907            1.64493406684822
1908            sage: R(-2).zeta()
1909            0.000000000000000
1910            sage: R(1).zeta()
1911            +infinity
1912
1913        Computing zeta using PARI is much more efficient in difficult cases.
1914        Here's how to compute zeta with at least a given precision:
1915       
1916             sage: z = pari.new_with_bits_prec(2, 53).zeta(); z
1917             1.644934066848226436472415167              # 32-bit
1918             1.6449340668482264364724151666460251892    # 64-bit
1919
1920        Note that the number of bits of precision in the constructor only
1921        effects the internel precision of the pari number, not the number
1922        of digits that gets displayed.  To increase that you must
1923        use \code{pari.set_real_precision}.
1924       
1925             sage: type(z)
1926             <type 'sage.libs.pari.gen.gen'>
1927             sage: R(z)
1928             1.64493406684822
1929        """
1930        cdef RealNumber x
1931        x = self._new()
1932        _sig_on
1933        mpfr_zeta(x.value, self.value, (<RealField>self._parent).rnd)
1934        _sig_off
1935        return x
1936
1937    def algdep(self, n):
1938        """
1939        Returns a polynomial of degree at most $n$ which is approximately
1940        satisfied by this number.  Note that the returned polynomial
1941        need not be irreducible, and indeed usually won't be if this number
1942        is a good approximation to an algebraic number of degree less than $n$.
1943
1944        ALGORITHM: Uses the PARI C-library algdep command.
1945
1946        EXAMPLE:
1947             sage: r = sqrt(2); r
1948             1.41421356237309
1949             sage: r.algdep(5)
1950             x^5 - x^4 - 2*x^3 + x^2 + 2      # 32-bit
1951             x^4 - 4*x^2 + 4                  # 64-bit
1952        """
1953        return sage.rings.arith.algdep(self,n)
1954
1955    def algebraic_dependency(self, n):
1956        """
1957         Returns a polynomial of degree at most $n$ which is approximately
1958         satisfied by this number.  Note that the returned polynomial
1959         need not be irreducible, and indeed usually won't be if this number
1960         is a good approximation to an algebraic number of degree less than $n$.
1961
1962         ALGORITHM: Uses the PARI C-library algdep command.
1963
1964         EXAMPLE:
1965              sage: r = sqrt(2); r
1966              1.41421356237309
1967              sage: r.algdep(5)
1968              x^5 - x^4 - 2*x^3 + x^2 + 2    # 32-bit
1969              x^4 - 4*x^2 + 4                # 64-bit
1970        """
1971        return sage.rings.arith.algdep(self,n)
1972
1973    def nth_root(self, int n):
1974        r"""
1975        Returns an $n^{th}$ root of self.
1976
1977        INPUT:
1978            n -- A positive number, rounded down to the nearest integer.
1979                 Note that $n$ should be less than $\code{sys.maxint}$.
1980       
1981        EXAMPLES:
1982            sage: R = RealField()
1983            sage: R(8).nth_root(3)
1984            2.00000000000000
1985            sage: R(8).nth_root(3.7)    # illustrate rounding down
1986            2.00000000000000
1987            sage: R(-8).nth_root(3)
1988            -2.00000000000000
1989            sage: R(0).nth_root(3)
1990            0.000000000000000
1991            sage: R(32).nth_root(-1)
1992            Traceback (most recent call last):
1993            ...
1994            ValueError: n must be nonnegative
1995            sage: R(32).nth_root(1.0)
1996            32.0000000000000
1997           
1998        Note that for negative numbers, any even root returns NaN
1999            sage: R(-2).nth_root(6)
2000            NaN
2001
2002        The $n^{th}$ root of 0 is defined to be 0, for any $n$
2003            sage: R(0).nth_root(6)
2004            0.000000000000000
2005
2006            sage: R(0).nth_root(7)
2007            0.000000000000000
2008       
2009        AUTHOR: Didier Deshommes (2007-02)
2010        REFEREE: David Harvey
2011
2012        TODO: (trac \#294) the underlying mpfr_root function is unforgivably
2013        slow when n is large. e.g. RealNumber(8).nth_root(100000) is very slow.
2014        This should be investigated further and possibly discussed with the
2015        mpfr developers.
2016        """
2017        cdef RealNumber x
2018       
2019        if n < 0:
2020            raise ValueError, "n must be nonnegative"
2021       
2022        x = self._new()
2023        _sig_on
2024        mpfr_root(x.value, self.value, n, (<RealField>self._parent).rnd)
2025        _sig_off
2026        return x
2027
2028   
2029RR = RealField()
2030
2031
2032def create_RealNumber(s, int base=10, int pad=0, rnd="RNDN", min_prec=53):
2033    r"""
2034    Return the real number defined by the string s as an element of
2035    \code{RealField(prec=n)}, where n potentially has slightly more
2036    (controlled by pad) bits than given by s.
2037
2038    INPUT:
2039        s -- a string that defines a real number (or something whose
2040             string representation defines a number)
2041        base -- an integer between 2 and 36
2042        pad -- an integer >= 0.
2043        rnd -- rounding mode: RNDN, RNDZ, RNDU, RNDD
2044        min_prec -- number will have at least this many bits of precision, no matter what.
2045
2046    EXAMPLES:
2047        sage: RealNumber('2.3')
2048        2.29999999999999
2049        sage: RealNumber(10)
2050        10.0000000000000
2051        sage: RealNumber('1.0000000000000000000000000000000000')
2052        1.0000000000000000000000000000000000   
2053    """
2054    if not isinstance(s, str):
2055        s = str(s)
2056    if base == 10:
2057        bits = int(3.32192*len(s))
2058    else:
2059        bits = int(math.log(base,2)*len(s))
2060    R = RealField(prec=max(bits+pad, min_prec), rnd=rnd)
2061    return RealNumber(R, s, base)
2062
2063
2064
2065def is_RealField(x):
2066    return PY_TYPE_CHECK(x, RealField)
2067
2068def is_RealNumber(x):
2069    return PY_TYPE_CHECK(x, RealNumber)
2070
2071def __create__RealField_version0(prec, sci_not, rnd):
2072    return RealField(prec, sci_not, rnd)
2073
2074def __create__RealNumber_version0(parent, x, base=10):
2075    return RealNumber(parent, x, base=base)
Note: See TracBrowser for help on using the repository browser.