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