-
-
Notifications
You must be signed in to change notification settings - Fork 405
/
real_arb.pyx
2705 lines (2247 loc) · 91.3 KB
/
real_arb.pyx
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
r"""
Arbitrary precision real intervals using Arb
AUTHORS:
- Clemens Heuberger (2014-10-21): Initial version.
This is a binding to the optional `Arb library <http://fredrikj.net/arb/>`_; it
may be useful to refer to its documentation for more details.
You may have to run ``sage -i arb`` to use the arb library.
Parts of the documentation for this module are copied or adapted from
Arb's own documentation, licenced under the GNU General Public License
version 2, or later.
Comparison
==========
.. WARNING::
Identical :class:`RealBall` objects are understood to give
permission for algebraic simplification. This assumption is made
to improve performance. For example, setting ``z = x*x`` sets `z`
to a ball enclosing the set `\{t^2 : t \in x\}` and not the
(generally larger) set `\{tu : t \in x, u \in x\}`.
Two elements are equal if and only if they are the same object
or if both are exact and equal::
sage: from sage.rings.real_arb import RealBallField # optional - arb
sage: RBF = RealBallField() # optional - arb
sage: a = RBF(1) # optional - arb
sage: b = RBF(1) # optional - arb
sage: a is b # optional - arb
False
sage: a == b # optional - arb
True
sage: a = RBF(1/3) # optional - arb
sage: b = RBF(1/3) # optional - arb
sage: a.is_exact() # optional - arb
False
sage: b.is_exact() # optional - arb
False
sage: a is b # optional - arb
False
sage: a == b # optional - arb
False
A ball is non-zero if and only if it does not contain zero. ::
sage: a = RBF(RIF(-0.5, 0.5)) # optional - arb
sage: bool(a) # optional - arb
False
sage: a != 0 # optional - arb
False
sage: b = RBF(1/3) # optional - arb
sage: bool(b) # optional - arb
True
sage: b != 0 # optional - arb
True
A ball ``left`` is less than a ball ``right`` if all elements of
``left`` are less than all elements of ``right``. ::
sage: a = RBF(RIF(1, 2)) # optional - arb
sage: b = RBF(RIF(3, 4)) # optional - arb
sage: a < b # optional - arb
True
sage: a <= b # optional - arb
True
sage: a > b # optional - arb
False
sage: a >= b # optional - arb
False
sage: a = RBF(RIF(1, 3)) # optional - arb
sage: b = RBF(RIF(2, 4)) # optional - arb
sage: a < b # optional - arb
False
sage: a <= b # optional - arb
False
sage: a > b # optional - arb
False
sage: a >= b # optional - arb
False
Comparisons with Sage symbolic infinities work with some limitations::
sage: -infinity < RBF(1) < +infinity # optional - arb
True
sage: -infinity < RBF(infinity) # optional - arb
True
sage: RBF(infinity) < infinity # optional - arb
False
sage: RBF(NaN) < infinity # optional - arb
Traceback (most recent call last):
...
ValueError: infinite but not with +/- phase
sage: 1/RBF(0) <= infinity # optional - arb
Traceback (most recent call last):
...
ValueError: infinite but not with +/- phase
Comparisons between elements of real ball fields, however, support special
values and should be preferred::
sage: RBF(NaN) < RBF(infinity) # optional - arb
False
sage: 1/RBF(0) <= RBF(infinity) # optional - arb
True
TESTS::
sage: from sage.rings.real_arb import RBF # optional - arb
sage: (RBF(pi) * identity_matrix(QQ, 3)).parent() # optional - arb
Full MatrixSpace of 3 by 3 dense matrices over Real ball field with 53 bits precision
Classes and Methods
===================
"""
#*****************************************************************************
# Copyright (C) 2014 Clemens Heuberger <clemens.heuberger@aau.at>
#
# Distributed under the terms of the GNU General Public License (GPL)
# as published by the Free Software Foundation; either version 2 of
# the License, or (at your option) any later version.
# http://www.gnu.org/licenses/
#*****************************************************************************
include 'sage/ext/interrupt.pxi'
include "sage/ext/python.pxi"
import operator
import sage.symbolic.constants
from sage.rings.integer_ring import ZZ
from sage.rings.rational_field import QQ
from sage.rings.real_mpfi import RealIntervalField, RealIntervalField_class
from sage.structure.unique_representation import UniqueRepresentation
cimport sage.rings.integer
cimport sage.rings.rational
cimport sage.structure.element
from sage.libs.arb.arb cimport *
from sage.libs.arb.arf cimport arf_t, arf_init, arf_get_mpfr, arf_set_mpfr, arf_clear, arf_set_mag, arf_set
from sage.libs.arb.arf cimport arf_equal, arf_is_nan, arf_is_neg_inf, arf_is_pos_inf, arf_get_mag
from sage.libs.arb.mag cimport mag_t, mag_init, mag_clear, mag_add, mag_set_d, MAG_BITS, mag_is_inf, mag_is_finite, mag_zero
from sage.libs.flint.flint cimport flint_free
from sage.libs.flint.fmpz cimport fmpz_t, fmpz_init, fmpz_get_mpz, fmpz_set_mpz, fmpz_clear
from sage.libs.flint.fmpq cimport fmpq_t, fmpq_init, fmpq_set_mpq, fmpq_clear
from sage.libs.gmp.mpz cimport mpz_fits_ulong_p, mpz_fits_slong_p, mpz_get_ui, mpz_get_si
from sage.libs.mpfi cimport mpfi_get_left, mpfi_get_right, mpfi_interv_fr
from sage.libs.mpfr cimport mpfr_t, mpfr_init2, mpfr_clear, mpfr_sgn, MPFR_PREC_MIN
from sage.libs.mpfr cimport GMP_RNDN, GMP_RNDU, GMP_RNDD, GMP_RNDZ
from sage.rings.real_double cimport RealDoubleElement
from sage.rings.real_mpfr cimport RealField_class, RealField, RealNumber
from sage.structure.element cimport Element, ModuleElement, RingElement
cdef void mpfi_to_arb(arb_t target, const mpfi_t source, const long precision):
"""
Convert an MPFI interval to an Arb ball.
INPUT:
- ``target`` -- an ``arb_t``.
- ``source`` -- an ``mpfi_t``.
- ``precision`` -- an integer `\ge 2`.
"""
cdef mpfr_t left
cdef mpfr_t right
if _do_sig(precision): sig_on()
mpfr_init2(left, precision)
mpfr_init2(right, precision)
mpfi_get_left(left, source)
mpfi_get_right(right, source)
arb_set_interval_mpfr(target,
left,
right,
precision)
mpfr_clear(left)
mpfr_clear(right)
cdef int arb_to_mpfi(mpfi_t target, arb_t source, const long precision) except -1:
"""
Convert an Arb ball to an MPFI interval.
INPUT:
- ``target`` -- an ``mpfi_t``.
- ``source`` -- an ``arb_t``.
- ``precision`` -- an integer `\ge 2`.
EXAMPLES::
sage: from sage.rings.real_arb import RBF # optional - arb
sage: RIF(RBF(2)**(2**100)) # optional - arb, indirect doctest
Traceback (most recent call last):
...
ArithmeticError: Error converting arb to mpfi. Overflow?
"""
cdef mpfr_t left
cdef mpfr_t right
mpfr_init2(left, precision)
mpfr_init2(right, precision)
try:
sig_on()
arb_get_interval_mpfr(left, right, source)
mpfi_interv_fr(target, left, right)
sig_off()
except RuntimeError:
raise ArithmeticError("Error converting arb to mpfi. Overflow?")
finally:
mpfr_clear(left)
mpfr_clear(right)
class RealBallField(UniqueRepresentation, Parent):
r"""
An approximation of the field of real numbers using mid-rad intervals, also
known as balls.
INPUT:
- ``precision`` -- an integer `\ge 2`.
EXAMPLES::
sage: from sage.rings.real_arb import RealBallField # optional - arb
sage: RBF = RealBallField() # optional - arb; indirect doctest
sage: RBF(1) # optional - arb
1.000000000000000
::
sage: from sage.rings.real_arb import RBF # optional - arb
sage: (1/2*RBF(1)) + AA(sqrt(2)) - 1 + polygen(QQ, x) # optional - arb
x + [0.914213562373095 +/- 4.10e-16]
TESTS::
sage: RBF.bracket(RBF(1/2), RBF(1/3)) # optional - arb
[+/- 5.56e-17]
sage: RBF.cardinality() # optional - arb
+Infinity
sage: RBF.cartesian_product(QQ).an_element()**2 # optional - arb
([1.440000000000000 +/- 4.98e-16], 1/4)
sage: RBF.coerce_embedding() is None # optional - arb
True
sage: loads(dumps(RBF)) is RBF # optional - arb
True
sage: RBF['x'].gens_dict_recursive() # optional - arb
{'x': x}
sage: RBF.is_finite() # optional - arb
False
sage: RBF.is_zero() # optional - arb
False
sage: RBF.one() # optional - arb
1.000000000000000
sage: RBF.zero() # optional - arb
0
"""
Element = RealBall
@staticmethod
def __classcall__(cls, long precision=53, category=None):
r"""
Normalize the arguments for caching.
TESTS::
sage: from sage.rings.real_arb import RealBallField # optional - arb
sage: RealBallField(53) is RealBallField() # optional - arb
True
"""
return super(RealBallField, cls).__classcall__(cls, precision, category)
def __init__(self, precision, category):
r"""
Initialize the real ball field.
INPUT:
- ``precision`` -- an integer `\ge 2`.
EXAMPLES::
sage: from sage.rings.real_arb import RealBallField # optional - arb
sage: RBF = RealBallField() # optional - arb
sage: RBF(1) # optional - arb
1.000000000000000
sage: RealBallField(0) # optional - arb
Traceback (most recent call last):
...
ValueError: Precision must be at least 2.
sage: RealBallField(1) # optional - arb
Traceback (most recent call last):
...
ValueError: Precision must be at least 2.
"""
if precision < 2:
raise ValueError("Precision must be at least 2.")
super(RealBallField, self).__init__(
#category=category or sage.categories.magmas_and_additive_magmas.MagmasAndAdditiveMagmas().Infinite(),
# FIXME: RBF is not even associative, but CompletionFunctor only works with rings.
category=category or sage.categories.rings.Rings().Infinite())
self._prec = precision
def _repr_(self):
r"""
String representation of ``self``.
EXAMPLES::
sage: from sage.rings.real_arb import RealBallField # optional - arb
sage: RealBallField() # optional - arb
Real ball field with 53 bits precision
sage: RealBallField(106) # optional - arb
Real ball field with 106 bits precision
"""
return "Real ball field with {} bits precision".format(self._prec)
def _coerce_map_from_(self, other):
r"""
Parents that canonically coerce into real ball fields include:
- some exact or lazy parents representing subsets of the reals, such as
``ZZ``, ``QQ``, ``AA``, and ``RLF``;
- real ball fields with a larger precision.
TESTS::
sage: from sage.rings.real_arb import RealBallField # optional - arb
sage: RealBallField().has_coerce_map_from(RealBallField(54)) # optional - arb
True
sage: RealBallField().has_coerce_map_from(RealBallField(52)) # optional - arb
False
sage: RealBallField().has_coerce_map_from(RIF) # optional - arb
False
sage: RealBallField().has_coerce_map_from(SR) # optional - arb
False
sage: RealBallField().has_coerce_map_from(RR) # optional - arb
False
"""
from sage.rings.qqbar import AA
from sage.rings.real_lazy import RLF
if isinstance(other, RealBallField):
return (other._prec >= self._prec)
elif (other is ZZ) or (other is QQ) or (other is AA) or (other is RLF):
return True
else:
return False
def _element_constructor_(self, mid=None, rad=None):
"""
Convert ``mid`` to an element of this real ball field, perhaps
non-canonically.
In addition to the inputs supported by
:meth:`ElementConstructor.__init__`,
anything that is convertible to a real interval can also be used to
construct a real ball::
sage: from sage.rings.real_arb import RBF # optional - arb
sage: RBF(RIF(0, 1)) # optional - arb; indirect doctest
[+/- 1.01]
sage: RBF(1) # optional - arb
1.000000000000000
sage: RBF(x) # optional - arb
Traceback (most recent call last):
...
TypeError: unable to convert x to a RealIntervalFieldElement
Various symbolic constants can be converted without going through real
intervals. (This is faster and yields tighter error bounds.) ::
sage: RBF(e) # optional - arb
[2.718281828459045 +/- 5.49e-16]
sage: RBF(pi) # optional - arb
[3.141592653589793 +/- 5.62e-16]
"""
try:
return self.element_class(self, mid, rad)
except TypeError:
pass
try:
return self.element_class(self, mid.pyobject(), rad)
except (AttributeError, TypeError):
pass
try:
mid = RealIntervalField(self._prec)(mid)
except TypeError:
raise TypeError("unable to convert {} to a RealIntervalFieldElement".format(mid))
return self.element_class(self, mid, rad)
def gens(self):
r"""
EXAMPLE::
sage: from sage.rings.real_arb import RBF # optional - arb
sage: RBF.gens() # optional - arb
(1.000000000000000,)
sage: RBF.gens_dict() # optional - arb
{'1.000000000000000': 1.000000000000000}
"""
return (self.one(),)
def base(self):
"""
Real ball fields are their own base.
EXAMPLE::
sage: from sage.rings.real_arb import RBF # optional - arb
sage: RBF.base() # optional - arb
Real ball field with 53 bits precision
"""
return self
def base_ring(self):
"""
Real ball fields are their own base ring.
EXAMPLE::
sage: from sage.rings.real_arb import RBF # optional - arb
sage: RBF.base_ring() # optional - arb
Real ball field with 53 bits precision
"""
return self
def construction(self):
"""
Return the construction of a real ball field as a completion of the
rationals.
EXAMPLES::
sage: from sage.rings.real_arb import RealBallField # optional - arb
sage: RBF = RealBallField(42) # optional - arb
sage: functor, base = RBF.construction() # optional - arb
sage: functor, base # optional - arb
(Completion[+Infinity], Rational Field)
sage: functor(base) is RBF # optional - arb
True
"""
from sage.categories.pushout import CompletionFunctor
functor = CompletionFunctor(sage.rings.infinity.Infinity,
self._prec,
{'type': 'Ball'})
return functor, QQ
def precision(self):
"""
Return the bit precision used for operations on elements of this field.
EXAMPLES::
sage: from sage.rings.real_arb import RealBallField # optional - arb
sage: RealBallField().precision() # optional - arb
53
"""
return self._prec
def is_exact(self):
"""
Real ball fields are not exact.
EXAMPLES::
sage: from sage.rings.real_arb import RealBallField # optional - arb
sage: RealBallField().is_exact() # optional - arb
False
"""
return False
def characteristic(self):
"""
Real ball fields have characteristic zero.
EXAMPLES::
sage: from sage.rings.real_arb import RealBallField # optional - arb
sage: RealBallField().characteristic() # optional - arb
0
"""
return 0
def some_elements(self):
"""
Real ball fields contain exact balls, inexact balls, infinities, and
more.
EXAMPLES::
sage: from sage.rings.real_arb import RBF # optional - arb
sage: RBF.some_elements() # optional - arb
[1.000000000000000,
[0.3333333333333333 +/- 7.04e-17],
[-4.733045976388941e+363922934236666733021124 +/- 3.46e+363922934236666733021108],
[+/- inf],
[+/- inf],
nan]
"""
return [self(1), self(1)/3,
-self(2)**(sage.rings.integer.Integer(2)**80),
self(sage.rings.infinity.Infinity), ~self(0),
self.element_class(self, sage.symbolic.constants.NotANumber())]
# Ball functions of non-ball arguments
def sinpi(self, x):
"""
Return a ball enclosing sin(πx).
This works even if ``x`` itself is not a ball, and may be faster or
more accurate where ``x`` is a rational number.
.. seealso:: :meth:`~sage.rings.real_arb.RealBall.sin`
EXAMPLES::
sage: from sage.rings.real_arb import RBF # optional - arb
sage: RBF.sinpi(1) # optional - arb
0
sage: RBF.sinpi(1/3) # optional - arb
[0.866025403784439 +/- 5.15e-16]
sage: RBF.sinpi(1 + 2^(-100)) # optional - arb
[-2.478279624546525e-30 +/- 5.92e-46]
TESTS::
sage: RBF.sinpi(RLF(sqrt(2))) # optional - arb
[-0.96390253284988 +/- 4.11e-15]
"""
cdef RealBall res, x_as_ball
cdef sage.rings.rational.Rational x_as_Rational
cdef fmpq_t tmpq
res = self.element_class(self)
try:
x_as_Rational = QQ.coerce(x)
try:
if _do_sig(self._prec): sig_on()
fmpq_init(tmpq)
fmpq_set_mpq(tmpq, x_as_Rational.value)
arb_sin_pi_fmpq(res.value, tmpq, self._prec)
if _do_sig(self._prec): sig_off()
finally:
fmpq_clear(tmpq)
return res
except TypeError:
pass
x_as_ball = self.coerce(x)
if _do_sig(self._prec): sig_on()
arb_sin_pi(res.value, x_as_ball.value, self._prec)
if _do_sig(self._prec): sig_off()
return res
def cospi(self, x):
"""
Return a ball enclosing cos(πx).
This works even if ``x`` itself is not a ball, and may be faster or
more accurate where ``x`` is a rational number.
.. seealso:: :meth:`~sage.rings.real_arb.RealBall.cos`
EXAMPLES::
sage: from sage.rings.real_arb import RBF # optional - arb
sage: RBF.cospi(1) # optional - arb
-1.000000000000000
sage: RBF.cospi(1/3) # optional - arb
0.5000000000000000
TESTS::
sage: RBF.cospi(RLF(sqrt(2))) # optional - arb
[-0.26625534204142 +/- 5.38e-15]
"""
cdef RealBall res, x_as_ball
cdef sage.rings.rational.Rational x_as_Rational
cdef fmpq_t tmpq
res = self.element_class(self)
try:
x_as_Rational = QQ.coerce(x)
try:
if _do_sig(self._prec): sig_on()
fmpq_init(tmpq)
fmpq_set_mpq(tmpq, x_as_Rational.value)
arb_cos_pi_fmpq(res.value, tmpq, self._prec)
if _do_sig(self._prec): sig_off()
finally:
fmpq_clear(tmpq)
return res
except TypeError:
pass
x_as_ball = self.coerce(x)
if _do_sig(self._prec): sig_on()
arb_cos_pi(res.value, x_as_ball.value, self._prec)
if _do_sig(self._prec): sig_off()
return res
def gamma(self, x):
"""
Return a ball enclosing the gamma function of ``x``.
This works even if ``x`` itself is not a ball, and may be more
efficient in the case where ``x`` is an integer or a rational number.
.. seealso:: :meth:`~sage.rings.real_arb.RealBall.gamma`
EXAMPLES::
sage: from sage.rings.real_arb import RBF # optional - arb
sage: RBF.gamma(5) # optional - arb
24.00000000000000
sage: RBF.gamma(10**20) # optional - arb
[+/- 5.92e+1956570551809674821757]
sage: RBF.gamma(1/3) # optional - arb
[2.678938534707747 +/- 9.06e-16]
sage: RBF.gamma(-5) # optional - arb
nan
TESTS::
sage: RBF.gamma(RLF(pi)) # optional - arb
[2.2880377953400 +/- 4.29e-14]
"""
cdef RealBall res
cdef sage.rings.integer.Integer x_as_Integer
cdef sage.rings.rational.Rational x_as_Rational
cdef fmpz_t tmpz
cdef fmpq_t tmpq
res = self.element_class(self)
try:
x_as_Integer = ZZ.coerce(x)
try:
if _do_sig(self._prec): sig_on()
fmpz_init(tmpz)
fmpz_set_mpz(tmpz, x_as_Integer.value)
arb_gamma_fmpz(res.value, tmpz, self._prec)
if _do_sig(self._prec): sig_off()
finally:
fmpz_clear(tmpz)
return res
except TypeError:
pass
try:
x_as_Rational = QQ.coerce(x)
try:
if _do_sig(self._prec): sig_on()
fmpq_init(tmpq)
fmpq_set_mpq(tmpq, x_as_Rational.value)
arb_gamma_fmpq(res.value, tmpq, self._prec)
if _do_sig(self._prec): sig_off()
finally:
fmpq_clear(tmpq)
return res
except TypeError:
pass
return self.coerce(x).gamma()
def zeta(self, s):
"""
Return a ball enclosing the Riemann zeta function of ``s``.
This works even if ``s`` itself is not a ball, and may be more
efficient in the case where ``s`` is an integer.
.. SEEALSO:: :meth:`~sage.rings.real_arb.RealBall.zeta`
EXAMPLES::
sage: from sage.rings.real_arb import RBF # optional - arb
sage: RBF.zeta(3) # abs tol 5e-16, optional - arb
[1.202056903159594 +/- 2.87e-16]
sage: RBF.zeta(1) # optional - arb
nan
sage: RBF.zeta(1/2) # optional - arb
[-1.460354508809587 +/- 1.94e-16]
"""
cdef RealBall res
cdef sage.rings.integer.Integer s_as_Integer
try:
s_as_Integer = ZZ.coerce(s)
if mpz_fits_ulong_p(s_as_Integer.value):
res = self.element_class(self)
if _do_sig(self._prec): sig_on()
arb_zeta_ui(res.value, mpz_get_ui(s_as_Integer.value), self._prec)
if _do_sig(self._prec): sig_off()
return res
except TypeError:
pass
return self.coerce(s).zeta()
def bernoulli(self, n):
"""
Return a ball enclosing the ``n``-th Bernoulli number.
EXAMPLES::
sage: from sage.rings.real_arb import RBF # optional - arb
sage: [RBF.bernoulli(n) for n in range(4)] # optional - arb
[1.000000000000000, -0.5000000000000000, [0.1666666666666667 +/- 7.04e-17], 0]
sage: RBF.bernoulli(2**20) # optional - arb
[-1.823002872104961e+5020717 +/- 7.16e+5020701]
sage: RBF.bernoulli(2**1000) # optional - arb
Traceback (most recent call last):
...
ValueError: argument too large
TESTS::
sage: RBF.bernoulli(2r) # optional - arb
[0.1666666666666667 +/- 7.04e-17]
sage: RBF.bernoulli(2/3) # optional - arb
Traceback (most recent call last):
...
TypeError: no canonical coercion from Rational Field to Integer Ring
sage: RBF.bernoulli(-1) # optional - arb
Traceback (most recent call last):
...
ValueError: expected a nonnegative index
"""
cdef RealBall res
cdef sage.rings.integer.Integer n_as_Integer = ZZ.coerce(n)
if mpz_fits_ulong_p(n_as_Integer.value):
res = self.element_class(self)
if _do_sig(self._prec): sig_on()
arb_bernoulli_ui(res.value, mpz_get_ui(n_as_Integer.value), self._prec)
if _do_sig(self._prec): sig_off()
return res
elif n_as_Integer < 0:
raise ValueError("expected a nonnegative index")
else:
# TODO: Fall back to a Sage implementation in this case?
raise ValueError("argument too large")
def fibonacci(self, n):
"""
Return a ball enclosing the ``n``-th Fibonacci number.
EXAMPLES::
sage: from sage.rings.real_arb import RBF # optional - arb
sage: [RBF.fibonacci(n) for n in xrange(7)] # optional - arb
[0,
1.000000000000000,
1.000000000000000,
2.000000000000000,
3.000000000000000,
5.000000000000000,
8.000000000000000]
sage: RBF.fibonacci(-2) # optional - arb
-1.000000000000000
sage: RBF.fibonacci(10**20) # optional - arb
[3.78202087472056e+20898764024997873376 +/- 4.01e+20898764024997873361]
"""
cdef fmpz_t tmpz
cdef RealBall res = self.element_class(self)
cdef sage.rings.integer.Integer n_as_Integer = ZZ.coerce(n)
try:
if _do_sig(self._prec): sig_on()
fmpz_init(tmpz)
fmpz_set_mpz(tmpz, n_as_Integer.value)
arb_fib_fmpz(res.value, tmpz, self._prec)
if _do_sig(self._prec): sig_off()
finally:
fmpz_clear(tmpz)
return res
cdef inline bint _do_sig(long prec):
"""
Whether signal handlers should be installed for calls to arb.
TESTS::
sage: from sage.rings.real_arb import RealBallField # optional - arb
sage: _ = RealBallField()(1).psi() # optional - arb; indirect doctest
sage: _ = RealBallField(1500)(1).psi() # optional - arb
"""
return (prec > 1000)
cdef inline long prec(RealBall ball):
return ball._parent._prec
cdef class RealBall(RingElement):
"""
Hold one ``arb_t`` of the `Arb library
<http://fredrikj.net/arb/>`_
EXAMPLES::
sage: from sage.rings.real_arb import RealBallField # optional - arb
sage: a = RealBallField()(RIF(1)) # optional - arb; indirect doctest
sage: b = a.psi() # optional - arb
sage: b # optional - arb
[-0.577215664901533 +/- 3.85e-16]
sage: RIF(b) # optional - arb
-0.577215664901533?
"""
def __cinit__(self):
"""
Allocate memory for the encapsulated value.
EXAMPLES::
sage: from sage.rings.real_arb import RealBallField # optional - arb
sage: RealBallField()(RIF(1)) # optional - arb; indirect doctest
1.000000000000000
"""
arb_init(self.value)
def __dealloc__(self):
"""
Deallocate memory of the encapsulated value.
EXAMPLES::
sage: from sage.rings.real_arb import RealBallField # optional - arb
sage: a = RealBallField()(RIF(1)) # optional - arb; indirect doctest
sage: del a # optional - arb
"""
arb_clear(self.value)
def __init__(self, parent, mid=None, rad=None):
"""
Initialize the :class:`RealBall`.
INPUT:
- ``parent`` -- a :class:`RealBallField`.
- ``mid`` (optional) -- ball midpoint, see examples below. If omitted,
initialize the ball to zero, ignoring the ``rad`` argument.
- ``rad`` (optional) -- a :class:`RealNumber` or a Python float, ball
radius. If the midpoint is not exactly representable in
floating-point, the radius is adjusted to account for the roundoff
error.
EXAMPLES::
sage: from sage.rings.real_arb import RealBallField # optional - arb
sage: RBF = RealBallField() # optional - arb
sage: RBF() # optional - arb
0
One can create exact real balls using elements of various exact parents,
or using floating-point numbers::
sage: RBF(3) # optional - arb
3.000000000000000
sage: RBF(3r) # optional - arb
3.000000000000000
sage: RBF(1/3) # optional - arb
[0.3333333333333333 +/- 7.04e-17]
sage: RBF(3.14) # optional - arb
[3.140000000000000 +/- 1.25e-16]
::
sage: RBF(3, 0.125) # optional - arb
[3e+0 +/- 0.126]
sage: RBF(pi, 0.125r) # optional - arb
[3e+0 +/- 0.267]
Note that integers and floating-point numbers are ''not'' rounded to
the parent's precision::
sage: b = RBF(11111111111111111111111111111111111111111111111); b # optional - arb
[1.111111111111111e+46 +/- 1.12e+30]
sage: b.mid().exact_rational() # optional - arb
11111111111111111111111111111111111111111111111
Similarly, converting a real ball from one real ball field to another
(with a different precision) only changes the way it is displayed and
the precision of operations involving it, not the actual representation
of its center::
sage: RBF100 = RealBallField(100) # optional - arb
sage: b100 = RBF100(1/3); b100 # optional - arb
[0.333333333333333333333333333333 +/- 4.65e-31]
sage: b53 = RBF(b100); b53 # optional - arb
[0.3333333333333333 +/- 3.34e-17]
sage: RBF100(b53) # optional - arb
[0.333333333333333333333333333333 +/- 4.65e-31]
Special values are supported::
sage: RBF(oo).mid(), RBF(-oo).mid(), RBF(unsigned_infinity).mid() # optional - arb
(+infinity, -infinity, 0.000000000000000)
sage: RBF(NaN) # optional - arb
nan
TESTS::
sage: from sage.rings.real_arb import RealBall # optional - arb
sage: RealBall(RBF, sage.symbolic.constants.Pi()) # abs tol 1e-16, optional - arb
[3.141592653589793 +/- 5.62e-16]
sage: RealBall(RBF, sage.symbolic.constants.Log2()) # abs tol 1e-16, optional - arb
[0.693147180559945 +/- 4.06e-16]
sage: RealBall(RBF, sage.symbolic.constants.Catalan()) # optional - arb
[0.915965594177219 +/- 9.43e-17]
sage: RealBall(RBF, sage.symbolic.constants.Khinchin()) # optional - arb
[2.685452001065306 +/- 6.82e-16]
sage: RealBall(RBF, sage.symbolic.constants.Glaisher()) # optional - arb
[1.282427129100623 +/- 6.02e-16]
sage: RealBall(RBF, sage.symbolic.constants.e) # optional - arb
[2.718281828459045 +/- 5.49e-16]
"""
cdef fmpz_t tmpz
cdef fmpq_t tmpq
cdef arf_t tmpr
cdef mag_t tmpm
Element.__init__(self, parent)
if mid is None:
return
elif isinstance(mid, RealBall):
arb_set(self.value, (<RealBall> mid).value) # no rounding!
elif isinstance(mid, int):
arb_set_si(self.value, PyInt_AS_LONG(mid)) # no rounding!
elif isinstance(mid, sage.rings.integer.Integer):
if _do_sig(prec(self)): sig_on()
fmpz_init(tmpz)
fmpz_set_mpz(tmpz, (<sage.rings.integer.Integer> mid).value)
arb_set_fmpz(self.value, tmpz) # no rounding!
fmpz_clear(tmpz)
if _do_sig(prec(self)): sig_off()
elif isinstance(mid, sage.rings.rational.Rational):
if _do_sig(prec(self)): sig_on()
fmpq_init(tmpq)
fmpq_set_mpq(tmpq, (<sage.rings.rational.Rational> mid).value)
arb_set_fmpq(self.value, tmpq, prec(self))
fmpq_clear(tmpq)
if _do_sig(prec(self)): sig_off()
elif isinstance(mid, RealNumber):
if _do_sig(prec(self)): sig_on()
arf_init(tmpr)
arf_set_mpfr(tmpr, (<RealNumber> mid).value)
arb_set_arf(self.value, tmpr) # no rounding!
arf_clear(tmpr)
if _do_sig(prec(self)): sig_off()
elif isinstance(mid, sage.rings.infinity.AnInfinity):
if isinstance(mid, sage.rings.infinity.PlusInfinity):
arb_pos_inf(self.value)
elif isinstance(mid, sage.rings.infinity.MinusInfinity):
arb_neg_inf(self.value)
else:
arb_zero_pm_inf(self.value)
elif isinstance(mid, sage.symbolic.constants.Constant):
if _do_sig(prec(self)): sig_on()
try:
if isinstance(mid, sage.symbolic.constants.NotANumber):
arb_indeterminate(self.value)
elif isinstance(mid, sage.symbolic.constants.Pi):
arb_const_pi(self.value, prec(self))
elif isinstance(mid, sage.symbolic.constants.Log2):
arb_const_log2(self.value, prec(self))
elif isinstance(mid, sage.symbolic.constants.Catalan):
arb_const_catalan(self.value, prec(self))
elif isinstance(mid, sage.symbolic.constants.Khinchin):
arb_const_khinchin(self.value, prec(self))
elif isinstance(mid, sage.symbolic.constants.Glaisher):
arb_const_glaisher(self.value, prec(self))
else:
raise TypeError("unsupported constant")
finally:
if _do_sig(prec(self)): sig_off()
elif isinstance(mid, sage.symbolic.constants_c.E):
if _do_sig(prec(self)): sig_on()
arb_const_e(self.value, prec(self))
if _do_sig(prec(self)): sig_off()
elif isinstance(mid, RealIntervalFieldElement):
mpfi_to_arb(self.value,
(<RealIntervalFieldElement> mid).value,
prec(self))
else:
raise TypeError("unsupported midpoint type")