| [1324] | 1 | """ |
|---|
| [2437] | 2 | Double Precision Real Numbers |
|---|
| [1324] | 3 | |
|---|
| [1636] | 4 | EXAMPLES: |
|---|
| 5 | |
|---|
| 6 | We create the real double vector space of dimension $3$: |
|---|
| 7 | sage: V = RDF^3; V |
|---|
| 8 | Vector space of dimension 3 over Real Double Field |
|---|
| 9 | |
|---|
| 10 | Notice that this space is unique. |
|---|
| 11 | sage: V is RDF^3 |
|---|
| 12 | True |
|---|
| 13 | sage: V is FreeModule(RDF, 3) |
|---|
| 14 | True |
|---|
| 15 | sage: V is VectorSpace(RDF, 3) |
|---|
| 16 | True |
|---|
| 17 | |
|---|
| 18 | Also, you can instantly create a space of large dimension. |
|---|
| 19 | sage: V = RDF^10000 |
|---|
| [1324] | 20 | """ |
|---|
| [1534] | 21 | |
|---|
| [1283] | 22 | include '../ext/cdefs.pxi' |
|---|
| [1820] | 23 | include '../ext/stdsage.pxi' |
|---|
| [1283] | 24 | include '../ext/interrupt.pxi' |
|---|
| 25 | include '../gsl/gsl.pxi' |
|---|
| 26 | |
|---|
| [2581] | 27 | import math, operator |
|---|
| [2664] | 28 | from random import random |
|---|
| [1283] | 29 | |
|---|
| 30 | from sage.misc.sage_eval import sage_eval |
|---|
| 31 | |
|---|
| [1749] | 32 | import sage.rings.complex_double |
|---|
| 33 | import sage.rings.complex_field |
|---|
| 34 | |
|---|
| 35 | import sage.rings.integer |
|---|
| 36 | import sage.rings.rational |
|---|
| 37 | |
|---|
| [1848] | 38 | cdef class RealDoubleField_class(Field): |
|---|
| [1283] | 39 | """ |
|---|
| 40 | The field of real double precision numbers. |
|---|
| 41 | |
|---|
| [1820] | 42 | EXAMPLES: |
|---|
| 43 | sage: RR == RDF |
|---|
| 44 | False |
|---|
| 45 | sage: RDF == RealDoubleField() # RDF is the shorthand |
|---|
| 46 | True |
|---|
| [1283] | 47 | """ |
|---|
| 48 | |
|---|
| [2123] | 49 | def is_exact(self): |
|---|
| 50 | return False |
|---|
| 51 | |
|---|
| [2581] | 52 | def _latex_(self): |
|---|
| 53 | return "\\R" |
|---|
| 54 | |
|---|
| [1283] | 55 | def __repr__(self): |
|---|
| 56 | """ |
|---|
| 57 | Print out this real double field. |
|---|
| 58 | |
|---|
| 59 | EXAMPLES: |
|---|
| [1820] | 60 | sage: RealDoubleField() |
|---|
| [1283] | 61 | Real Double Field |
|---|
| 62 | sage: RDF |
|---|
| 63 | Real Double Field |
|---|
| 64 | """ |
|---|
| 65 | return "Real Double Field" |
|---|
| [1598] | 66 | |
|---|
| [2618] | 67 | def __cmp__(self, x): |
|---|
| 68 | """ |
|---|
| 69 | EXAMPLES: |
|---|
| 70 | sage: RDF == 5 |
|---|
| 71 | False |
|---|
| 72 | sage: loads(dumps(RDF)) == RDF |
|---|
| 73 | True |
|---|
| 74 | """ |
|---|
| 75 | if PY_TYPE_CHECK(x, RealDoubleField_class): |
|---|
| 76 | return 0 |
|---|
| 77 | return cmp(type(self), type(x)) |
|---|
| 78 | |
|---|
| [1283] | 79 | def __call__(self, x): |
|---|
| 80 | """ |
|---|
| 81 | Create a real double using x. |
|---|
| 82 | |
|---|
| 83 | EXAMPLES: |
|---|
| 84 | sage: RDF(1) |
|---|
| 85 | 1.0 |
|---|
| 86 | sage: RDF(2/3) |
|---|
| 87 | 0.666666666667 |
|---|
| 88 | |
|---|
| 89 | A TypeError is raised if the coercion doesn't make sense: |
|---|
| 90 | sage: RDF(QQ['x'].0) |
|---|
| 91 | Traceback (most recent call last): |
|---|
| 92 | ... |
|---|
| 93 | TypeError: cannot coerce nonconstant polynomial to float |
|---|
| 94 | |
|---|
| 95 | One can convert back and forth between double precision real |
|---|
| 96 | numbers and higher-precision ones, though of course there may |
|---|
| 97 | be loss of precision: |
|---|
| 98 | sage: a = RealField(200)(2).sqrt(); a |
|---|
| [1749] | 99 | 1.4142135623730950488016887242096980785696718753769480731766 |
|---|
| [1283] | 100 | sage: b = RDF(a); b |
|---|
| [1319] | 101 | 1.41421356237 |
|---|
| [1283] | 102 | sage: a.parent()(b) |
|---|
| [2552] | 103 | 1.4142135623700000000000000000000000000000000000000000000000 |
|---|
| [1283] | 104 | """ |
|---|
| 105 | return RealDoubleElement(x) |
|---|
| 106 | |
|---|
| [1848] | 107 | cdef _coerce_c_impl(self, x): |
|---|
| [1749] | 108 | """ |
|---|
| 109 | Canonical coercion of x to the real double field. |
|---|
| 110 | |
|---|
| 111 | The rings that canonically coerce to the real double field are: |
|---|
| 112 | * the real double field itself |
|---|
| 113 | * int, long, integer, and rational rings |
|---|
| 114 | * real mathematical constants |
|---|
| 115 | * the mpfr real field |
|---|
| 116 | |
|---|
| 117 | EXAMPLES: |
|---|
| 118 | sage: RDF._coerce_(5) |
|---|
| 119 | 5.0 |
|---|
| 120 | sage: RDF._coerce_(9499294r) |
|---|
| 121 | 9499294.0 |
|---|
| 122 | sage: RDF._coerce_(61/3) |
|---|
| 123 | 20.3333333333 |
|---|
| 124 | sage: parent(RDF(3) + CDF(5)) |
|---|
| 125 | Complex Double Field |
|---|
| 126 | sage: parent(CDF(5) + RDF(3)) |
|---|
| 127 | Complex Double Field |
|---|
| 128 | """ |
|---|
| 129 | if isinstance(x, (int, long, sage.rings.integer.Integer, |
|---|
| 130 | sage.rings.rational.Rational)): |
|---|
| 131 | return self(x) |
|---|
| 132 | import real_mpfr |
|---|
| 133 | return self._coerce_try(x, [sage.functions.constants.ConstantRing, |
|---|
| 134 | real_mpfr.RR]) |
|---|
| 135 | |
|---|
| 136 | |
|---|
| [1283] | 137 | def gen(self, n=0): |
|---|
| 138 | """ |
|---|
| 139 | Return the generator of the real double field. |
|---|
| 140 | EXAMPLES: |
|---|
| [1319] | 141 | sage: RDF.0 |
|---|
| [1283] | 142 | 1.0 |
|---|
| [1319] | 143 | sage: RDF.gens() |
|---|
| [1283] | 144 | (1.0,) |
|---|
| 145 | """ |
|---|
| 146 | if n != 0: |
|---|
| 147 | raise ValueError, "only 1 generator" |
|---|
| [1319] | 148 | return RealDoubleElement(1) |
|---|
| [1283] | 149 | |
|---|
| 150 | def ngens(self): |
|---|
| 151 | return 1 |
|---|
| 152 | |
|---|
| 153 | def is_atomic_repr(self): |
|---|
| 154 | """ |
|---|
| 155 | Returns True, to signify that elements of this field print |
|---|
| 156 | without sums, so parenthesis aren't required, e.g., in |
|---|
| 157 | coefficients of polynomials. |
|---|
| 158 | |
|---|
| 159 | EXAMPLES: |
|---|
| [1851] | 160 | sage: RDF.is_atomic_repr() |
|---|
| [1283] | 161 | True |
|---|
| 162 | """ |
|---|
| 163 | return True |
|---|
| 164 | |
|---|
| 165 | def is_finite(self): |
|---|
| 166 | """ |
|---|
| 167 | Returns False, since the field of real numbers is not finite. |
|---|
| 168 | Technical note: There exists an upper bound on the double representation. |
|---|
| 169 | |
|---|
| 170 | EXAMPLES: |
|---|
| [1851] | 171 | sage: RDF.is_finite() |
|---|
| [1283] | 172 | False |
|---|
| 173 | """ |
|---|
| 174 | return False |
|---|
| 175 | |
|---|
| 176 | def characteristic(self): |
|---|
| 177 | """ |
|---|
| 178 | Returns 0, since the field of real numbers has characteristic 0. |
|---|
| 179 | |
|---|
| 180 | EXAMPLES: |
|---|
| [1851] | 181 | sage: RDF.characteristic() |
|---|
| [1283] | 182 | 0 |
|---|
| 183 | """ |
|---|
| 184 | return 0 |
|---|
| [2664] | 185 | |
|---|
| 186 | cdef _new_c(self, double value): |
|---|
| 187 | cdef RealDoubleElement x |
|---|
| 188 | x = PY_NEW(RealDoubleElement) |
|---|
| 189 | x._value = value |
|---|
| 190 | return x |
|---|
| 191 | |
|---|
| 192 | def random_element(self, float min=-1, float max=1): |
|---|
| 193 | """ |
|---|
| 194 | Return a random element of this real double field in the interval [min, max]. |
|---|
| 195 | |
|---|
| 196 | EXAMPLES: |
|---|
| 197 | """ |
|---|
| 198 | return self._new_c((max-min)*random() + min) |
|---|
| [1283] | 199 | |
|---|
| 200 | def name(self): |
|---|
| 201 | return "RealDoubleField" |
|---|
| 202 | |
|---|
| 203 | def __hash__(self): |
|---|
| [1319] | 204 | return 1455926870 #return hash(self.name()) |
|---|
| [1283] | 205 | |
|---|
| 206 | def pi(self): |
|---|
| 207 | """ |
|---|
| 208 | Returns pi to double-precision. |
|---|
| 209 | |
|---|
| 210 | EXAMPLES: |
|---|
| 211 | sage: R = RealField(100) |
|---|
| 212 | sage: RDF.pi() |
|---|
| 213 | 3.14159265359 |
|---|
| 214 | sage: RDF.pi().sqrt()/2 |
|---|
| [1319] | 215 | 0.886226925453 |
|---|
| [1283] | 216 | """ |
|---|
| 217 | return self(M_PI) |
|---|
| 218 | |
|---|
| 219 | |
|---|
| 220 | def euler_constant(self): |
|---|
| 221 | """ |
|---|
| 222 | Returns Euler's gamma constant to double precision |
|---|
| 223 | |
|---|
| 224 | EXAMPLES: |
|---|
| 225 | sage: RDF.euler_constant() |
|---|
| 226 | 0.577215664902 |
|---|
| 227 | """ |
|---|
| 228 | return self(M_EULER) |
|---|
| 229 | |
|---|
| 230 | def log2(self): |
|---|
| 231 | """ |
|---|
| 232 | Returns log(2) to the precision of this field. |
|---|
| 233 | |
|---|
| 234 | EXAMPLES: |
|---|
| 235 | sage: RDF.log2() |
|---|
| 236 | 0.69314718056 |
|---|
| 237 | sage: RDF(2).log() |
|---|
| 238 | 0.69314718056 |
|---|
| 239 | """ |
|---|
| 240 | return self(M_LN2) |
|---|
| 241 | |
|---|
| 242 | def factorial(self, int n): |
|---|
| 243 | """ |
|---|
| 244 | Return the factorial of the integer n as a real number. |
|---|
| 245 | EXAMPLES: |
|---|
| 246 | sage: RDF.factorial(100) |
|---|
| 247 | 9.33262154439e+157 |
|---|
| 248 | """ |
|---|
| 249 | if n < 0: |
|---|
| 250 | raise ArithmeticError, "n must be nonnegative" |
|---|
| 251 | return self(gsl_sf_fact(n)) |
|---|
| 252 | |
|---|
| 253 | def zeta(self, n=2): |
|---|
| 254 | """ |
|---|
| 255 | Return an $n$-th root of unity in the real field, |
|---|
| 256 | if one exists, or raise a ValueError otherwise. |
|---|
| 257 | |
|---|
| 258 | EXAMPLES: |
|---|
| 259 | sage: RDF.zeta() |
|---|
| 260 | -1.0 |
|---|
| 261 | sage: RDF.zeta(1) |
|---|
| 262 | 1.0 |
|---|
| [1319] | 263 | sage: RDF.zeta(5) |
|---|
| [1283] | 264 | Traceback (most recent call last): |
|---|
| 265 | ... |
|---|
| 266 | ValueError: No 5th root of unity in self |
|---|
| 267 | """ |
|---|
| 268 | if n == 1: |
|---|
| 269 | return self(1) |
|---|
| 270 | elif n == 2: |
|---|
| 271 | return self(-1) |
|---|
| 272 | raise ValueError, "No %sth root of unity in self"%n |
|---|
| 273 | |
|---|
| 274 | |
|---|
| 275 | |
|---|
| [2321] | 276 | def new_RealDoubleElement(): |
|---|
| 277 | cdef RealDoubleElement x |
|---|
| 278 | x = PY_NEW(RealDoubleElement) |
|---|
| 279 | return x |
|---|
| [1283] | 280 | |
|---|
| [1820] | 281 | cdef class RealDoubleElement(FieldElement): |
|---|
| [2664] | 282 | def __new__(self, x=None): |
|---|
| 283 | (<Element>self)._parent = _RDF |
|---|
| 284 | |
|---|
| [1283] | 285 | def __init__(self, x): |
|---|
| 286 | self._value = float(x) |
|---|
| [1820] | 287 | |
|---|
| [2618] | 288 | def __reduce__(self): |
|---|
| 289 | """ |
|---|
| 290 | EXAMPLES: |
|---|
| 291 | sage: a = RDF(-2.7) |
|---|
| 292 | sage: loads(dumps(a)) == a |
|---|
| 293 | True |
|---|
| 294 | """ |
|---|
| 295 | return RealDoubleElement, (self._value, ) |
|---|
| 296 | |
|---|
| [1820] | 297 | cdef _new_c(self, double value): |
|---|
| 298 | cdef RealDoubleElement x |
|---|
| 299 | x = PY_NEW(RealDoubleElement) |
|---|
| 300 | x._value = value |
|---|
| 301 | return x |
|---|
| [1283] | 302 | |
|---|
| 303 | def real(self): |
|---|
| 304 | """ |
|---|
| 305 | Returns itself -- we're already real. |
|---|
| 306 | |
|---|
| 307 | EXAMPLES: |
|---|
| 308 | sage: a = RDF(3) |
|---|
| 309 | sage: a.real() |
|---|
| 310 | 3.0 |
|---|
| 311 | """ |
|---|
| 312 | return self |
|---|
| 313 | |
|---|
| 314 | def imag(self): |
|---|
| 315 | """ |
|---|
| 316 | Returns the imaginary part of this number. (hint: it's zero.) |
|---|
| 317 | |
|---|
| 318 | EXAMPLES: |
|---|
| 319 | sage: a = RDF(3) |
|---|
| 320 | sage: a.imag() |
|---|
| 321 | 0.0 |
|---|
| 322 | """ |
|---|
| [1319] | 323 | return RealDoubleElement(0) |
|---|
| [1283] | 324 | |
|---|
| 325 | def __complex__(self): |
|---|
| 326 | """ |
|---|
| 327 | EXAMPLES: |
|---|
| 328 | sage: a = 2303 |
|---|
| 329 | sage: RDF(a) |
|---|
| 330 | 2303.0 |
|---|
| 331 | sage: complex(RDF(a)) |
|---|
| 332 | (2303+0j) |
|---|
| 333 | """ |
|---|
| 334 | return complex(self._value,0) |
|---|
| 335 | |
|---|
| 336 | def parent(self): |
|---|
| 337 | """ |
|---|
| 338 | Return the real double field, which is the parent of self. |
|---|
| 339 | |
|---|
| 340 | EXAMPLES: |
|---|
| [1319] | 341 | sage: a = RDF(2.3) |
|---|
| [1283] | 342 | sage: a.parent() |
|---|
| 343 | Real Double Field |
|---|
| 344 | sage: parent(a) |
|---|
| 345 | Real Double Field |
|---|
| 346 | """ |
|---|
| 347 | return RDF |
|---|
| 348 | |
|---|
| 349 | def __repr__(self): |
|---|
| 350 | """ |
|---|
| 351 | Return print version of self. |
|---|
| 352 | |
|---|
| 353 | EXAMPLES: |
|---|
| 354 | sage: a = RDF(2); a |
|---|
| 355 | 2.0 |
|---|
| 356 | sage: a^2 |
|---|
| 357 | 4.0 |
|---|
| 358 | """ |
|---|
| 359 | return self.str() |
|---|
| 360 | |
|---|
| [2581] | 361 | def _latex_(self): # todo -- this is terrible if sci not. |
|---|
| [1283] | 362 | return self.str() |
|---|
| 363 | |
|---|
| 364 | def __hash__(self): |
|---|
| [2618] | 365 | return 1455926870 |
|---|
| 366 | #return hash(self.str()) |
|---|
| [1283] | 367 | |
|---|
| 368 | def _im_gens_(self, codomain, im_gens): |
|---|
| 369 | return codomain(self) # since 1 |--> 1 |
|---|
| 370 | |
|---|
| [2581] | 371 | def str(self): |
|---|
| 372 | """ |
|---|
| 373 | Return string representation of self. |
|---|
| 374 | |
|---|
| 375 | EXAMPLES: |
|---|
| 376 | sage: a = RDF('4.5'); a.str() |
|---|
| 377 | '4.5' |
|---|
| 378 | sage: a = RDF('49203480923840.2923904823048'); a.str() |
|---|
| 379 | '4.92034809238e+13' |
|---|
| 380 | sage: a = RDF(1)/RDF(0); a.str() |
|---|
| 381 | 'inf' |
|---|
| 382 | sage: a = -RDF(1)/RDF(0); a.str() |
|---|
| 383 | '-inf' |
|---|
| 384 | sage: a = RDF(0)/RDF(0); a.str() |
|---|
| 385 | 'nan' |
|---|
| 386 | """ |
|---|
| [1283] | 387 | if gsl_isnan(self._value): |
|---|
| 388 | return "nan" |
|---|
| 389 | else: |
|---|
| 390 | v = gsl_isinf(self._value) |
|---|
| 391 | if v == 1: |
|---|
| 392 | return "inf" |
|---|
| 393 | elif v == -1: |
|---|
| 394 | return "-inf" |
|---|
| [2581] | 395 | return str(self._value) |
|---|
| [1283] | 396 | |
|---|
| [2581] | 397 | def __copy__(self): |
|---|
| 398 | """ |
|---|
| 399 | Return copy of self, which since self is immutable, is just self. |
|---|
| 400 | |
|---|
| 401 | EXAMPLES: |
|---|
| 402 | sage: r = RDF('-1.6') |
|---|
| 403 | sage: r.__copy__() is r |
|---|
| 404 | True |
|---|
| 405 | """ |
|---|
| 406 | return self |
|---|
| 407 | #cdef RealDoubleElement z |
|---|
| 408 | #z = RealDoubleElement(self._value) |
|---|
| 409 | #return z |
|---|
| [1283] | 410 | |
|---|
| 411 | def integer_part(self): |
|---|
| 412 | """ |
|---|
| 413 | If in decimal this number is written n.defg, returns n. |
|---|
| [2581] | 414 | |
|---|
| 415 | EXAMPLES: |
|---|
| 416 | sage: r = RDF('-1.6') |
|---|
| 417 | sage: a = r.integer_part(); a |
|---|
| 418 | -1 |
|---|
| 419 | sage: type(a) |
|---|
| 420 | <type 'sage.rings.integer.Integer'> |
|---|
| [1283] | 421 | """ |
|---|
| [2581] | 422 | return sage.rings.integer.Integer(int(self._value)) |
|---|
| [1283] | 423 | |
|---|
| 424 | |
|---|
| 425 | ######################## |
|---|
| 426 | # Basic Arithmetic |
|---|
| 427 | ######################## |
|---|
| 428 | def __invert__(self): |
|---|
| 429 | """ |
|---|
| 430 | Compute the multiplicative inverse of self. |
|---|
| 431 | |
|---|
| 432 | EXAMPLES: |
|---|
| 433 | sage: a = RDF(-1.5)*RDF(2.5) |
|---|
| 434 | sage: a.__invert__() |
|---|
| 435 | -0.266666666667 |
|---|
| 436 | """ |
|---|
| 437 | return RealDoubleElement(1/self._value) |
|---|
| 438 | |
|---|
| [1820] | 439 | cdef ModuleElement _add_c_impl(self, ModuleElement right): |
|---|
| 440 | """ |
|---|
| 441 | Add two real numbers with the same parent. |
|---|
| 442 | |
|---|
| 443 | EXAMPLES: |
|---|
| [1851] | 444 | sage: RDF('-1.5') + RDF('2.5') |
|---|
| [1820] | 445 | 1.0 |
|---|
| 446 | """ |
|---|
| 447 | return self._new_c(self._value + (<RealDoubleElement>right)._value) |
|---|
| 448 | |
|---|
| 449 | cdef ModuleElement _sub_c_impl(self, ModuleElement right): |
|---|
| [1283] | 450 | """ |
|---|
| 451 | Subtract two real numbers with the same parent. |
|---|
| 452 | |
|---|
| 453 | EXAMPLES: |
|---|
| [1851] | 454 | sage: RDF('-1.5') - RDF('2.5') |
|---|
| [1283] | 455 | -4.0 |
|---|
| 456 | """ |
|---|
| [1820] | 457 | return self._new_c(self._value - (<RealDoubleElement>right)._value) |
|---|
| [1283] | 458 | |
|---|
| [1820] | 459 | cdef RingElement _mul_c_impl(self, RingElement right): |
|---|
| [1283] | 460 | """ |
|---|
| 461 | Multiply two real numbers with the same parent. |
|---|
| 462 | |
|---|
| 463 | EXAMPLES: |
|---|
| [1851] | 464 | sage: RDF('-1.5') * RDF('2.5') |
|---|
| [1283] | 465 | -3.75 |
|---|
| 466 | """ |
|---|
| [1820] | 467 | return self._new_c(self._value * (<RealDoubleElement>right)._value) |
|---|
| [1283] | 468 | |
|---|
| [1820] | 469 | cdef RingElement _div_c_impl(self, RingElement right): |
|---|
| [1283] | 470 | """ |
|---|
| 471 | EXAMPLES: |
|---|
| [1851] | 472 | sage: RDF('-1.5') / RDF('2.5') |
|---|
| [1283] | 473 | -0.6 |
|---|
| [1851] | 474 | sage: RDF(1)/RDF(0) |
|---|
| 475 | inf |
|---|
| [1283] | 476 | """ |
|---|
| [1820] | 477 | return self._new_c(self._value / (<RealDoubleElement>right)._value) |
|---|
| [1283] | 478 | |
|---|
| [2247] | 479 | cdef ModuleElement _neg_c_impl(self): |
|---|
| 480 | """ |
|---|
| 481 | Negates a real number. |
|---|
| 482 | |
|---|
| 483 | EXAMPLES: |
|---|
| 484 | sage: -RDF('-1.5') |
|---|
| 485 | 1.5 |
|---|
| 486 | """ |
|---|
| [1820] | 487 | return self._new_c(-self._value) |
|---|
| [1283] | 488 | |
|---|
| 489 | def __abs__(self): |
|---|
| 490 | return self.abs() |
|---|
| 491 | |
|---|
| 492 | cdef RealDoubleElement abs(RealDoubleElement self): |
|---|
| 493 | return RealDoubleElement(abs(self._value)) |
|---|
| 494 | |
|---|
| 495 | def __lshift__(x, y): |
|---|
| 496 | """ |
|---|
| 497 | LShifting a double is not supported; nor is lshifting a RealDoubleElement. |
|---|
| 498 | """ |
|---|
| 499 | raise TypeError, "unsupported operand type(s) for <<: '%s' and '%s'"%(typeof(self), typeof(n)) |
|---|
| 500 | |
|---|
| 501 | def __rshift__(x, y): |
|---|
| 502 | """ |
|---|
| 503 | RShifting a double is not supported; nor is rshifting a RealDoubleElement. |
|---|
| 504 | """ |
|---|
| 505 | raise TypeError, "unsupported operand type(s) for >>: '%s' and '%s'"%(typeof(self), typeof(n)) |
|---|
| 506 | |
|---|
| 507 | def multiplicative_order(self): |
|---|
| 508 | if self == 1: |
|---|
| 509 | return 1 |
|---|
| 510 | elif self == -1: |
|---|
| 511 | return -1 |
|---|
| 512 | return sage.rings.infinity.infinity |
|---|
| 513 | |
|---|
| 514 | def sign(self): |
|---|
| 515 | """ |
|---|
| 516 | Returns -1,0, or 1 if self is negative, zero, or positive; respectively. |
|---|
| 517 | |
|---|
| 518 | Examples: |
|---|
| 519 | sage: RDF(-1.5).sign() |
|---|
| 520 | -1 |
|---|
| 521 | sage: RDF(0).sign() |
|---|
| 522 | 0 |
|---|
| 523 | sage: RDF(2.5).sign() |
|---|
| 524 | 1 |
|---|
| 525 | """ |
|---|
| 526 | if not self._value: |
|---|
| 527 | return 0 |
|---|
| 528 | if self._value > 0: |
|---|
| 529 | return 1 |
|---|
| 530 | return -1 |
|---|
| 531 | |
|---|
| 532 | |
|---|
| 533 | ################### |
|---|
| 534 | # Rounding etc |
|---|
| 535 | ################### |
|---|
| 536 | |
|---|
| 537 | def round(self): |
|---|
| 538 | """ |
|---|
| 539 | Given real number x, rounds up if fractional part is greater than .5, |
|---|
| 540 | rounds down if fractional part is lesser than .5. |
|---|
| 541 | EXAMPLES: |
|---|
| 542 | sage: RDF(0.49).round() |
|---|
| 543 | 0.0 |
|---|
| 544 | sage: RDF(0.51).round() |
|---|
| 545 | 1.0 |
|---|
| 546 | """ |
|---|
| 547 | return RealDoubleElement(round(self._value)) |
|---|
| 548 | |
|---|
| 549 | def floor(self): |
|---|
| 550 | """ |
|---|
| 551 | Returns the floor of this number |
|---|
| 552 | |
|---|
| 553 | EXAMPLES: |
|---|
| 554 | sage: RDF(2.99).floor() |
|---|
| 555 | 2 |
|---|
| 556 | sage: RDF(2.00).floor() |
|---|
| 557 | 2 |
|---|
| 558 | sage: RDF(-5/2).floor() |
|---|
| 559 | -3 |
|---|
| 560 | """ |
|---|
| [2581] | 561 | return sage.rings.integer.Integer(int(math.floor(self._value))) |
|---|
| [1283] | 562 | |
|---|
| 563 | def ceil(self): |
|---|
| 564 | """ |
|---|
| 565 | Returns the ceiling of this number |
|---|
| 566 | |
|---|
| 567 | OUTPUT: |
|---|
| 568 | integer |
|---|
| 569 | |
|---|
| 570 | EXAMPLES: |
|---|
| 571 | sage: RDF(2.99).ceil() |
|---|
| [1319] | 572 | 3 |
|---|
| [1283] | 573 | sage: RDF(2.00).ceil() |
|---|
| [1319] | 574 | 2 |
|---|
| [1283] | 575 | sage: RDF(-5/2).ceil() |
|---|
| [1319] | 576 | -2 |
|---|
| [1283] | 577 | """ |
|---|
| [2581] | 578 | return sage.rings.integer.Integer(int(math.ceil(self._value))) |
|---|
| [1283] | 579 | |
|---|
| 580 | def ceiling(self): |
|---|
| 581 | return self.ceil() |
|---|
| 582 | |
|---|
| 583 | def trunc(self): |
|---|
| 584 | """ |
|---|
| 585 | Truncates this number (returns integer part). |
|---|
| 586 | |
|---|
| 587 | EXAMPLES: |
|---|
| 588 | sage: RDF(2.99).trunc() |
|---|
| 589 | 2.0 |
|---|
| 590 | sage: RDF(-2.00).trunc() |
|---|
| 591 | -2.0 |
|---|
| 592 | sage: RDF(0.00).trunc() |
|---|
| 593 | 0.0 |
|---|
| 594 | """ |
|---|
| 595 | return RealDoubleElement(int(self._value)) |
|---|
| 596 | |
|---|
| 597 | def frac(self): |
|---|
| 598 | """ |
|---|
| 599 | frac returns a real number > -1 and < 1. it satisfies the |
|---|
| 600 | relation: |
|---|
| 601 | x = x.trunc() + x.frac() |
|---|
| 602 | |
|---|
| 603 | EXAMPLES: |
|---|
| 604 | sage: RDF(2.99).frac() |
|---|
| 605 | 0.99 |
|---|
| 606 | sage: RDF(2.50).frac() |
|---|
| 607 | 0.5 |
|---|
| 608 | sage: RDF(-2.79).frac() |
|---|
| 609 | -0.79 |
|---|
| 610 | """ |
|---|
| [1820] | 611 | return self._new_c(self._value - int(self._value)) |
|---|
| [1283] | 612 | |
|---|
| 613 | ########################################### |
|---|
| 614 | # Conversions |
|---|
| 615 | ########################################### |
|---|
| 616 | |
|---|
| 617 | def __float__(self): |
|---|
| 618 | return self._value |
|---|
| 619 | |
|---|
| 620 | def __int__(self): |
|---|
| 621 | """ |
|---|
| 622 | Returns integer truncation of this real number. |
|---|
| 623 | """ |
|---|
| 624 | return int(self._value) |
|---|
| 625 | |
|---|
| 626 | def __long__(self): |
|---|
| 627 | """ |
|---|
| 628 | Returns long integer truncation of this real number. |
|---|
| 629 | """ |
|---|
| 630 | return long(self._value) |
|---|
| 631 | |
|---|
| 632 | def _complex_number_(self): |
|---|
| 633 | return sage.rings.complex_field.ComplexField()(self) |
|---|
| 634 | |
|---|
| 635 | def _complex_double_(self): |
|---|
| [1749] | 636 | return sage.rings.complex_double.ComplexDoubleField(self) |
|---|
| [1283] | 637 | |
|---|
| 638 | def _pari_(self): |
|---|
| 639 | return sage.libs.pari.all.pari.new_with_bits_prec("%.15e"%self._value, 64) |
|---|
| 640 | |
|---|
| 641 | |
|---|
| 642 | ########################################### |
|---|
| 643 | # Comparisons: ==, !=, <, <=, >, >= |
|---|
| 644 | ########################################### |
|---|
| 645 | |
|---|
| 646 | def is_NaN(self): |
|---|
| 647 | return bool(gsl_isnan(self._value)) |
|---|
| 648 | |
|---|
| [2581] | 649 | def is_positive_infinity(self): |
|---|
| 650 | return bool(gsl_isinf(self._value) > 0) |
|---|
| 651 | |
|---|
| 652 | def is_negative_infinity(self): |
|---|
| 653 | return bool(gsl_isinf(self._value) < 0) |
|---|
| 654 | |
|---|
| [2790] | 655 | def is_infinity(self): |
|---|
| 656 | return bool(gsl_isinf(self._value)) |
|---|
| 657 | |
|---|
| [1820] | 658 | def __richcmp__(left, right, int op): |
|---|
| 659 | return (<Element>left)._richcmp(right, op) |
|---|
| 660 | |
|---|
| 661 | cdef int _cmp_c_impl(left, Element right) except -2: |
|---|
| 662 | if left._value < (<RealDoubleElement>right)._value: |
|---|
| [1319] | 663 | return -1 |
|---|
| [1820] | 664 | elif left._value > (<RealDoubleElement>right)._value: |
|---|
| [1319] | 665 | return 1 |
|---|
| 666 | return 0 |
|---|
| [1283] | 667 | |
|---|
| 668 | |
|---|
| 669 | ############################ |
|---|
| 670 | # Special Functions |
|---|
| 671 | ############################ |
|---|
| 672 | |
|---|
| 673 | def sqrt(self): |
|---|
| 674 | """ |
|---|
| 675 | Return a square root of self. |
|---|
| 676 | |
|---|
| 677 | If self is negative a complex number is returned. |
|---|
| 678 | |
|---|
| 679 | If you use self.square_root() then a real number will always |
|---|
| 680 | be returned (though it will be NaN if self is negative). |
|---|
| 681 | |
|---|
| 682 | EXAMPLES: |
|---|
| 683 | sage: r = 4.0 |
|---|
| 684 | sage: r.sqrt() |
|---|
| [1749] | 685 | 2.00000000000000 |
|---|
| [1283] | 686 | sage: r.sqrt()^2 == r |
|---|
| 687 | True |
|---|
| 688 | |
|---|
| 689 | sage: r = 4344 |
|---|
| 690 | sage: r.sqrt() |
|---|
| [1749] | 691 | 65.9090282131363 |
|---|
| 692 | sage: r.sqrt()^2 - r |
|---|
| [2552] | 693 | 0.000000000000000 |
|---|
| [1283] | 694 | |
|---|
| 695 | sage: r = -2.0 |
|---|
| 696 | sage: r.sqrt() |
|---|
| [1749] | 697 | 1.41421356237309*I |
|---|
| [1283] | 698 | """ |
|---|
| 699 | if self >= 0: |
|---|
| 700 | return self.square_root() |
|---|
| 701 | return self._complex_double_().sqrt() |
|---|
| 702 | |
|---|
| 703 | |
|---|
| 704 | def square_root(self): |
|---|
| 705 | """ |
|---|
| 706 | Return a square root of self. A real number will always be |
|---|
| 707 | returned (though it will be NaN if self is negative). |
|---|
| 708 | |
|---|
| 709 | Use self.sqrt() to get a complex number if self is negative. |
|---|
| 710 | |
|---|
| 711 | EXAMPLES: |
|---|
| 712 | sage: r = -2.0 |
|---|
| 713 | sage: r.square_root() |
|---|
| 714 | NaN |
|---|
| 715 | sage: r.sqrt() |
|---|
| [1749] | 716 | 1.41421356237309*I |
|---|
| [1283] | 717 | """ |
|---|
| [1820] | 718 | return self._new_c(sqrt(self._value)) |
|---|
| [1283] | 719 | |
|---|
| 720 | def cube_root(self): |
|---|
| 721 | """ |
|---|
| 722 | Return the cubic root (defined over the real numbers) of self. |
|---|
| 723 | |
|---|
| 724 | EXAMPLES: |
|---|
| 725 | sage: r = RDF(125.0); r.cube_root() |
|---|
| 726 | 5.0 |
|---|
| 727 | sage: r = RDF(-119.0) |
|---|
| [1377] | 728 | sage: r.cube_root()^3 - r # output is random, depending on arch. |
|---|
| 729 | 0.0 |
|---|
| [1283] | 730 | """ |
|---|
| 731 | return self.nth_root(3) |
|---|
| 732 | |
|---|
| 733 | |
|---|
| 734 | def nth_root(self, int n): |
|---|
| 735 | """ |
|---|
| [1324] | 736 | Returns the $n^{th}$ root of self. |
|---|
| [1283] | 737 | EXAMPLES: |
|---|
| 738 | sage: r = RDF(-125.0); r.nth_root(3) |
|---|
| 739 | -5.0 |
|---|
| 740 | sage: r.nth_root(5) |
|---|
| 741 | -2.6265278044 |
|---|
| 742 | """ |
|---|
| 743 | if n == 0: |
|---|
| 744 | return RealDoubleElement(float('nan')) |
|---|
| 745 | if self._value < 0 and GSL_IS_EVEN(n): |
|---|
| 746 | pass #return self._complex_double_().pow(1.0/n) |
|---|
| 747 | else: |
|---|
| 748 | return RealDoubleElement(self.__nth_root(n)) |
|---|
| 749 | |
|---|
| 750 | cdef double __nth_root(RealDoubleElement self, int n): |
|---|
| 751 | cdef int m |
|---|
| 752 | cdef double x |
|---|
| 753 | cdef double x0 |
|---|
| 754 | cdef double dx |
|---|
| 755 | cdef double dx0 |
|---|
| 756 | m = n-1 |
|---|
| 757 | x = ( m + self._value ) / n |
|---|
| 758 | x0 = 0 |
|---|
| 759 | dx = abs(x - x0) |
|---|
| 760 | dx0= dx + 1 |
|---|
| 761 | while dx < dx0: |
|---|
| 762 | x0= x |
|---|
| 763 | dx0 = dx |
|---|
| 764 | x = ( m*x + self._value / gsl_pow_int(x,m) ) / n |
|---|
| 765 | dx=abs(x - x0) |
|---|
| 766 | return x |
|---|
| 767 | |
|---|
| 768 | def __pow(self, RealDoubleElement exponent): |
|---|
| [1820] | 769 | return self._new_c(self._value**exponent._value) |
|---|
| [1283] | 770 | |
|---|
| 771 | def __pow_int(self, int exponent): |
|---|
| [1820] | 772 | return self._new_c(gsl_pow_int(self._value, exponent)) |
|---|
| [1283] | 773 | |
|---|
| 774 | def __pow__(self, exponent, modulus): |
|---|
| 775 | """ |
|---|
| 776 | Compute self raised to the power of exponent, rounded in |
|---|
| 777 | the direction specified by the parent of self. |
|---|
| 778 | |
|---|
| 779 | If the result is not a real number, self and the exponent are |
|---|
| 780 | both coerced to complex numbers (with sufficient precision), |
|---|
| 781 | then the exponentiation is computed in the complex numbers. |
|---|
| 782 | Thus this function can return either a real or complex number. |
|---|
| 783 | |
|---|
| 784 | EXAMPLES: |
|---|
| 785 | sage: a = RDF('1.23456') |
|---|
| 786 | sage: a^20 |
|---|
| [1319] | 787 | 67.6462977039 |
|---|
| [1283] | 788 | sage: a^a |
|---|
| [1319] | 789 | 1.29711148178 |
|---|
| [1283] | 790 | """ |
|---|
| 791 | cdef RealDoubleElement x |
|---|
| 792 | if isinstance(self, RealDoubleElement): |
|---|
| 793 | return self.__pow(RealDoubleElement(exponent)) |
|---|
| 794 | if isinstance(exponent, (int,Integer)): |
|---|
| 795 | return self.__pow_int(int(exponent)) |
|---|
| 796 | elif not isinstance(exponent, RealDoubleElement): |
|---|
| 797 | x = RealDoubleElement(exponent) |
|---|
| 798 | else: |
|---|
| 799 | x = exponent |
|---|
| 800 | return self.__pow(x) |
|---|
| 801 | |
|---|
| 802 | |
|---|
| [2581] | 803 | def _log_base(self, double log_of_base): |
|---|
| [1283] | 804 | if self._value < 2: |
|---|
| 805 | if self._value == 0: |
|---|
| 806 | return -1./0 |
|---|
| 807 | if self._value < 0: |
|---|
| 808 | return 0./0 |
|---|
| [1820] | 809 | return self._new_c(gsl_sf_log_1plusx(self._value - 1) / log_of_base) |
|---|
| 810 | return self._new_c(gsl_sf_log(self._value) / log_of_base) |
|---|
| [1283] | 811 | |
|---|
| [2581] | 812 | def log(self, base=None): |
|---|
| [1283] | 813 | """ |
|---|
| 814 | EXAMPLES: |
|---|
| 815 | sage: RDF(2).log() |
|---|
| 816 | 0.69314718056 |
|---|
| 817 | sage: RDF(2).log(2) |
|---|
| 818 | 1.0 |
|---|
| 819 | sage: RDF(2).log(pi) |
|---|
| 820 | 0.605511561398 |
|---|
| 821 | sage: RDF(2).log(10) |
|---|
| 822 | 0.301029995664 |
|---|
| 823 | sage: RDF(2).log(1.5) |
|---|
| 824 | 1.70951129135 |
|---|
| 825 | sage: RDF(0).log() |
|---|
| 826 | -inf |
|---|
| 827 | sage: RDF(-1).log() |
|---|
| 828 | nan |
|---|
| 829 | """ |
|---|
| [2581] | 830 | if base is None: |
|---|
| 831 | return self._log_base(1) |
|---|
| [1283] | 832 | else: |
|---|
| [1319] | 833 | if isinstance(base, RealDoubleElement): |
|---|
| [2581] | 834 | return self._log_base(base._log_base(1)) |
|---|
| [1283] | 835 | else: |
|---|
| [2581] | 836 | return self._log_base(gsl_sf_log(float(base))) |
|---|
| [1283] | 837 | |
|---|
| 838 | def log2(self): |
|---|
| 839 | """ |
|---|
| 840 | Returns log to the base 2 of self |
|---|
| 841 | |
|---|
| 842 | EXAMPLES: |
|---|
| 843 | sage: r = RDF(16.0) |
|---|
| 844 | sage: r.log2() |
|---|
| 845 | 4.0 |
|---|
| 846 | |
|---|
| 847 | sage: r = RDF(31.9); r.log2() |
|---|
| 848 | 4.99548451888 |
|---|
| 849 | |
|---|
| 850 | """ |
|---|
| [1820] | 851 | return self._new_c(gsl_sf_log(self._value) / M_LN2) |
|---|
| [1283] | 852 | |
|---|
| 853 | |
|---|
| 854 | def log10(self): |
|---|
| 855 | """ |
|---|
| 856 | Returns log to the base 10 of self |
|---|
| 857 | |
|---|
| 858 | EXAMPLES: |
|---|
| [1319] | 859 | sage: r = RDF('16.0'); r.log10() |
|---|
| [1283] | 860 | 1.20411998266 |
|---|
| 861 | sage: r.log() / log(10) |
|---|
| 862 | 1.20411998266 |
|---|
| [1319] | 863 | sage: r = RDF('39.9'); r.log10() |
|---|
| [1283] | 864 | 1.60097289569 |
|---|
| 865 | """ |
|---|
| [1820] | 866 | return self._new_c(gsl_sf_log(self._value) / M_LN10) |
|---|
| [1283] | 867 | |
|---|
| 868 | def logpi(self): |
|---|
| 869 | """ |
|---|
| 870 | Returns log to the base pi of self |
|---|
| 871 | |
|---|
| 872 | EXAMPLES: |
|---|
| [1319] | 873 | sage: r = RDF(16); r.logpi() |
|---|
| 874 | 2.42204624559 |
|---|
| [1283] | 875 | sage: r.log() / log(pi) |
|---|
| [1319] | 876 | 2.42204624559 |
|---|
| 877 | sage: r = RDF('39.9'); r.logpi() |
|---|
| 878 | 3.22030233461 |
|---|
| [1283] | 879 | """ |
|---|
| [1820] | 880 | return self._new_c(gsl_sf_log(self._value) / M_LNPI) |
|---|
| [1283] | 881 | |
|---|
| 882 | def exp(self): |
|---|
| 883 | r""" |
|---|
| 884 | Returns $e^\code{self}$ |
|---|
| 885 | |
|---|
| 886 | EXAMPLES: |
|---|
| 887 | sage: r = RDF(0.0) |
|---|
| 888 | sage: r.exp() |
|---|
| 889 | 1.0 |
|---|
| 890 | |
|---|
| [1319] | 891 | sage: r = RDF('32.3') |
|---|
| [1283] | 892 | sage: a = r.exp(); a |
|---|
| 893 | 1.06588847275e+14 |
|---|
| 894 | sage: a.log() |
|---|
| 895 | 32.3 |
|---|
| 896 | |
|---|
| [1319] | 897 | sage: r = RDF('-32.3') |
|---|
| [1283] | 898 | sage: r.exp() |
|---|
| 899 | 9.3818445885e-15 |
|---|
| 900 | """ |
|---|
| [1820] | 901 | return self._new_c(gsl_sf_exp(self._value)) |
|---|
| [1283] | 902 | |
|---|
| 903 | def exp2(self): |
|---|
| 904 | """ |
|---|
| 905 | Returns $2^\code{self}$ |
|---|
| 906 | |
|---|
| 907 | EXAMPLES: |
|---|
| 908 | sage: r = RDF(0.0) |
|---|
| 909 | sage: r.exp2() |
|---|
| 910 | 1.0 |
|---|
| 911 | |
|---|
| 912 | sage: r = RDF(32.0) |
|---|
| 913 | sage: r.exp2() |
|---|
| 914 | 4294967296.0 |
|---|
| 915 | |
|---|
| 916 | sage: r = RDF(-32.3) |
|---|
| 917 | sage: r.exp2() |
|---|
| 918 | 1.89117248253e-10 |
|---|
| 919 | |
|---|
| 920 | """ |
|---|
| [1820] | 921 | return self._new_c(gsl_sf_exp(self._value * M_LN2)) |
|---|
| [1283] | 922 | |
|---|
| 923 | def exp10(self): |
|---|
| 924 | r""" |
|---|
| 925 | Returns $10^\code{self}$ |
|---|
| 926 | |
|---|
| 927 | EXAMPLES: |
|---|
| 928 | sage: r = RDF(0.0) |
|---|
| 929 | sage: r.exp10() |
|---|
| 930 | 1.0 |
|---|
| 931 | |
|---|
| 932 | sage: r = RDF(32.0) |
|---|
| 933 | sage: r.exp10() |
|---|
| 934 | 1e+32 |
|---|
| 935 | |
|---|
| 936 | sage: r = RDF(-32.3) |
|---|
| 937 | sage: r.exp10() |
|---|
| 938 | 5.01187233627e-33 |
|---|
| 939 | """ |
|---|
| [1820] | 940 | return self._new_c(gsl_sf_exp(self._value * M_LN10)) |
|---|
| [1283] | 941 | |
|---|
| 942 | def cos(self): |
|---|
| 943 | """ |
|---|
| 944 | Returns the cosine of this number |
|---|
| 945 | |
|---|
| 946 | EXAMPLES: |
|---|
| 947 | sage: t=RDF.pi()/2 |
|---|
| 948 | sage: t.cos() |
|---|
| 949 | 6.12323399574e-17 |
|---|
| 950 | """ |
|---|
| [1820] | 951 | return self._new_c(gsl_sf_cos(self._value)) |
|---|
| [1283] | 952 | |
|---|
| 953 | def sin(self): |
|---|
| 954 | """ |
|---|
| 955 | Returns the sine of this number |
|---|
| 956 | |
|---|
| 957 | EXAMPLES: |
|---|
| 958 | sage: RDF(2).sin() |
|---|
| 959 | 0.909297426826 |
|---|
| 960 | """ |
|---|
| [1820] | 961 | return self._new_c(gsl_sf_sin(self._value)) |
|---|
| [1283] | 962 | |
|---|
| 963 | def tan(self): |
|---|
| 964 | """ |
|---|
| 965 | Returns the tangent of this number |
|---|
| 966 | |
|---|
| 967 | EXAMPLES: |
|---|
| 968 | sage: q = RDF.pi()/3 |
|---|
| 969 | sage: q.tan() |
|---|
| 970 | 1.73205080757 |
|---|
| 971 | sage: q = RDF.pi()/6 |
|---|
| 972 | sage: q.tan() |
|---|
| 973 | 0.57735026919 |
|---|
| 974 | """ |
|---|
| [1820] | 975 | return self._new_c(tan(self._value)) |
|---|
| [1283] | 976 | |
|---|
| 977 | def sincos(self): |
|---|
| 978 | """ |
|---|
| 979 | Returns a pair consisting of the sine and cosine. |
|---|
| 980 | |
|---|
| 981 | EXAMPLES: |
|---|
| 982 | sage: t = RDF.pi()/6 |
|---|
| 983 | sage: t.sincos() |
|---|
| 984 | (0.5, 0.866025403784) |
|---|
| 985 | """ |
|---|
| 986 | return self.sin(), self.cos() |
|---|
| 987 | |
|---|
| 988 | def hypot(self, other): |
|---|
| [1820] | 989 | return self._new_c(gsl_sf_hypot(self._value, float(other))) |
|---|
| [1283] | 990 | |
|---|
| 991 | def acos(self): |
|---|
| 992 | """ |
|---|
| 993 | Returns the inverse cosine of this number |
|---|
| 994 | |
|---|
| 995 | EXAMPLES: |
|---|
| 996 | sage: q = RDF.pi()/3 |
|---|
| 997 | sage: i = q.cos() |
|---|
| 998 | sage: i.acos() == q |
|---|
| 999 | True |
|---|
| 1000 | """ |
|---|
| [1820] | 1001 | return self._new_c(acos(self._value)) |
|---|
| [1283] | 1002 | |
|---|
| 1003 | def asin(self): |
|---|
| 1004 | """ |
|---|
| 1005 | Returns the inverse sine of this number |
|---|
| 1006 | |
|---|
| 1007 | EXAMPLES: |
|---|
| 1008 | sage: q = RDF.pi()/5 |
|---|
| 1009 | sage: i = q.sin() |
|---|
| 1010 | sage: i.asin() == q |
|---|
| 1011 | True |
|---|
| 1012 | """ |
|---|
| [1820] | 1013 | return self._new_c(asin(self._value)) |
|---|
| [1283] | 1014 | |
|---|
| 1015 | def atan(self): |
|---|
| 1016 | """ |
|---|
| 1017 | Returns the inverse tangent of this number |
|---|
| 1018 | |
|---|
| 1019 | EXAMPLES: |
|---|
| 1020 | sage: q = RDF.pi()/5 |
|---|
| 1021 | sage: i = q.tan() |
|---|
| 1022 | sage: i.atan() == q |
|---|
| 1023 | True |
|---|
| 1024 | """ |
|---|
| [1820] | 1025 | return self._new_c(atan(self._value)) |
|---|
| [1283] | 1026 | |
|---|
| 1027 | |
|---|
| 1028 | def cosh(self): |
|---|
| 1029 | """ |
|---|
| 1030 | Returns the hyperbolic cosine of this number |
|---|
| 1031 | |
|---|
| 1032 | EXAMPLES: |
|---|
| 1033 | sage: q = RDF.pi()/12 |
|---|
| 1034 | sage: q.cosh() |
|---|
| 1035 | 1.0344656401 |
|---|
| 1036 | """ |
|---|
| [1820] | 1037 | return self._new_c(cosh(self._value)) |
|---|
| [1283] | 1038 | |
|---|
| 1039 | def sinh(self): |
|---|
| 1040 | """ |
|---|
| 1041 | Returns the hyperbolic sine of this number |
|---|
| 1042 | |
|---|
| 1043 | EXAMPLES: |
|---|
| 1044 | sage: q = RDF.pi()/12 |
|---|
| 1045 | sage: q.sinh() |
|---|
| 1046 | 0.264800227602 |
|---|
| 1047 | |
|---|
| 1048 | """ |
|---|
| [1820] | 1049 | return self._new_c(sinh(self._value)) |
|---|
| [1283] | 1050 | |
|---|
| 1051 | def tanh(self): |
|---|
| 1052 | """ |
|---|
| 1053 | Returns the hyperbolic tangent of this number |
|---|
| 1054 | |
|---|
| 1055 | EXAMPLES: |
|---|
| 1056 | sage: q = RDF.pi()/12 |
|---|
| 1057 | sage: q.tanh() |
|---|
| 1058 | 0.255977789246 |
|---|
| 1059 | """ |
|---|
| [1820] | 1060 | return self._new_c(tanh(self._value)) |
|---|
| [1283] | 1061 | |
|---|
| 1062 | def acosh(self): |
|---|
| 1063 | """ |
|---|
| 1064 | Returns the hyperbolic inverse cosine of this number |
|---|
| 1065 | |
|---|
| 1066 | EXAMPLES: |
|---|
| 1067 | sage: q = RDF.pi()/2 |
|---|
| 1068 | sage: i = q.cosh() ; i |
|---|
| 1069 | 2.50917847866 |
|---|
| 1070 | sage: i.acosh() == q |
|---|
| 1071 | True |
|---|
| 1072 | """ |
|---|
| [1820] | 1073 | return self._new_c(gsl_acosh(self._value)) |
|---|
| [1283] | 1074 | |
|---|
| 1075 | def asinh(self): |
|---|
| 1076 | """ |
|---|
| 1077 | Returns the hyperbolic inverse sine of this number |
|---|
| 1078 | |
|---|
| 1079 | EXAMPLES: |
|---|
| 1080 | sage: q = RDF.pi()/2 |
|---|
| 1081 | sage: i = q.sinh() ; i |
|---|
| 1082 | 2.30129890231 |
|---|
| 1083 | sage: i.asinh() == q |
|---|
| 1084 | True |
|---|
| 1085 | """ |
|---|
| [1820] | 1086 | return self._new_c(gsl_asinh(self._value)) |
|---|
| [1283] | 1087 | |
|---|
| 1088 | def atanh(self): |
|---|
| 1089 | """ |
|---|
| 1090 | Returns the hyperbolic inverse tangent of this number |
|---|
| 1091 | |
|---|
| 1092 | EXAMPLES: |
|---|
| 1093 | sage: q = RDF.pi()/2 |
|---|
| 1094 | sage: i = q.tanh() ; i |
|---|
| 1095 | 0.917152335667 |
|---|
| [1377] | 1096 | sage: i.atanh() - q # output is random, depending on arch. |
|---|
| [1319] | 1097 | -4.4408920985e-16 |
|---|
| [1283] | 1098 | """ |
|---|
| [1820] | 1099 | return self._new_c(gsl_atanh(self._value)) |
|---|
| [1283] | 1100 | |
|---|
| 1101 | def agm(self, other): |
|---|
| 1102 | """ |
|---|
| 1103 | Return the arithmetic-geometric mean of self and other. The |
|---|
| 1104 | arithmetic-geometric mean is the common limit of the sequences |
|---|
| 1105 | $u_n$ and $v_n$, where $u_0$ is self, $v_0$ is other, |
|---|
| 1106 | $u_{n+1}$ is the arithmetic mean of $u_n$ and $v_n$, and |
|---|
| 1107 | $v_{n+1}$ is the geometric mean of u_n and v_n. If any operand |
|---|
| 1108 | is negative, the return value is \code{NaN}. |
|---|
| 1109 | """ |
|---|
| 1110 | return RealDoubleElement(sage.rings.all.RR(self).agm(sage.rings.all.RR(other))) |
|---|
| 1111 | |
|---|
| 1112 | def erf(self): |
|---|
| 1113 | """ |
|---|
| 1114 | Returns the value of the error function on self. |
|---|
| 1115 | |
|---|
| 1116 | EXAMPLES: |
|---|
| 1117 | sage: RDF(6).erf() |
|---|
| 1118 | 1.0 |
|---|
| 1119 | """ |
|---|
| [1820] | 1120 | return self._new_c(gsl_sf_erf(self._value)) |
|---|
| [1283] | 1121 | |
|---|
| 1122 | def gamma(self): |
|---|
| 1123 | """ |
|---|
| 1124 | The Euler gamma function. Return gamma of self. |
|---|
| 1125 | |
|---|
| 1126 | EXAMPLES: |
|---|
| 1127 | sage: RDF(6).gamma() |
|---|
| 1128 | 120.0 |
|---|
| 1129 | sage: RDF(1.5).gamma() |
|---|
| 1130 | 0.886226925453 |
|---|
| 1131 | """ |
|---|
| [1820] | 1132 | return self._new_c(gsl_sf_gamma(self._value)) |
|---|
| [1283] | 1133 | |
|---|
| 1134 | def zeta(self): |
|---|
| 1135 | r""" |
|---|
| 1136 | Return the Riemann zeta function evaluated at this real number. |
|---|
| 1137 | |
|---|
| 1138 | \note{PARI is vastly more efficient at computing the Riemann zeta |
|---|
| 1139 | function. See the example below for how to use it.} |
|---|
| 1140 | |
|---|
| 1141 | EXAMPLES: |
|---|
| 1142 | sage: RDF(2).zeta() |
|---|
| 1143 | 1.64493406685 |
|---|
| 1144 | sage: RDF.pi()^2/6 |
|---|
| 1145 | 1.64493406685 |
|---|
| [1378] | 1146 | sage: RDF(-2).zeta() # slightly random-ish arch dependent output |
|---|
| [1319] | 1147 | -2.37378795339e-18 |
|---|
| [1283] | 1148 | sage: RDF(1).zeta() |
|---|
| 1149 | inf |
|---|
| 1150 | """ |
|---|
| [1319] | 1151 | if self._value == 1: |
|---|
| [1820] | 1152 | return self._new_c(1)/self._new_c(0) |
|---|
| 1153 | return self._new_c(gsl_sf_zeta(self._value)) |
|---|
| [1283] | 1154 | |
|---|
| 1155 | def algdep(self, n): |
|---|
| 1156 | """ |
|---|
| 1157 | Returns a polynomial of degree at most $n$ which is approximately |
|---|
| 1158 | satisfied by this number. Note that the returned polynomial |
|---|
| 1159 | need not be irreducible, and indeed usually won't be if this number |
|---|
| 1160 | is a good approximation to an algebraic number of degree less than $n$. |
|---|
| 1161 | |
|---|
| 1162 | ALGORITHM: Uses the PARI C-library algdep command. |
|---|
| 1163 | |
|---|
| 1164 | EXAMPLE: |
|---|
| [1319] | 1165 | sage: r = RDF(2).sqrt(); r |
|---|
| 1166 | 1.41421356237 |
|---|
| 1167 | sage: r.algdep(5) |
|---|
| 1168 | x^4 - 2*x^2 |
|---|
| [1283] | 1169 | """ |
|---|
| 1170 | return sage.rings.arith.algdep(self,n) |
|---|
| 1171 | |
|---|
| 1172 | def algebraic_dependency(self, n): |
|---|
| 1173 | """ |
|---|
| [1319] | 1174 | Returns a polynomial of degree at most $n$ which is approximately |
|---|
| 1175 | satisfied by this number. Note that the returned polynomial |
|---|
| 1176 | need not be irreducible, and indeed usually won't be if this number |
|---|
| 1177 | is a good approximation to an algebraic number of degree less than $n$. |
|---|
| [1283] | 1178 | |
|---|
| [1319] | 1179 | ALGORITHM: Uses the PARI C-library algdep command. |
|---|
| [1283] | 1180 | |
|---|
| [1319] | 1181 | EXAMPLE: |
|---|
| 1182 | sage: r = sqrt(RDF(2)); r |
|---|
| 1183 | 1.41421356237 |
|---|
| 1184 | sage: r.algdep(5) |
|---|
| 1185 | x^4 - 2*x^2 |
|---|
| [1283] | 1186 | """ |
|---|
| 1187 | return sage.rings.arith.algdep(self,n) |
|---|
| 1188 | |
|---|
| 1189 | |
|---|
| 1190 | ##################################################### |
|---|
| 1191 | # unique objects |
|---|
| 1192 | ##################################################### |
|---|
| [1820] | 1193 | cdef RealDoubleField_class _RDF |
|---|
| 1194 | _RDF = RealDoubleField_class() |
|---|
| [1283] | 1195 | |
|---|
| [1820] | 1196 | RDF = _RDF # external interface |
|---|
| [1283] | 1197 | |
|---|
| [1820] | 1198 | def RealDoubleField(): |
|---|
| 1199 | global _RDF |
|---|
| 1200 | return _RDF |
|---|
| [1283] | 1201 | |
|---|
| [2618] | 1202 | |
|---|
| [2633] | 1203 | def is_RealDoubleElement(x): |
|---|
| 1204 | return PY_TYPE_CHECK(x, RealDoubleElement) |
|---|
| 1205 | |
|---|