Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagesmc
Path: blob/master/src/sage/rings/complex_field.py
8817 views
1
r"""
2
Field of Arbitrary Precision Complex Numbers
3
4
AUTHORS:
5
6
- William Stein (2006-01-26): complete rewrite
7
8
- Niles Johnson (2010-08): :trac:`3893`: ``random_element()`` should pass on
9
``*args`` and ``**kwds``.
10
11
- Travis Scrimshaw (2012-10-18): Added documentation for full coverage.
12
"""
13
14
#################################################################################
15
# Copyright (C) 2006 William Stein <[email protected]>
16
#
17
# Distributed under the terms of the GNU General Public License (GPL)
18
#
19
# http://www.gnu.org/licenses/
20
#*****************************************************************************
21
22
import complex_number
23
import complex_double
24
import field
25
import integer
26
import real_mpfr
27
import weakref
28
from sage.misc.sage_eval import sage_eval
29
30
from sage.structure.parent import Parent
31
from sage.structure.parent_gens import ParentWithGens
32
33
NumberField_quadratic = None
34
NumberFieldElement_quadratic = None
35
AlgebraicNumber_base = None
36
AlgebraicNumber = None
37
AlgebraicReal = None
38
AA = None
39
QQbar = None
40
SR = None
41
CDF = CLF = RLF = None
42
43
def late_import():
44
"""
45
Import the objects/modules after build (when needed).
46
47
TESTS::
48
49
sage: sage.rings.complex_field.late_import()
50
"""
51
global NumberField_quadratic
52
global NumberFieldElement_quadratic
53
global AlgebraicNumber_base
54
global AlgebraicNumber
55
global AlgebraicReal
56
global AA, QQbar, SR
57
global CLF, RLF, CDF
58
if NumberFieldElement_quadratic is None:
59
import sage.rings.number_field.number_field
60
import sage.rings.number_field.number_field_element_quadratic as nfeq
61
NumberField_quadratic = sage.rings.number_field.number_field.NumberField_quadratic
62
NumberFieldElement_quadratic = nfeq.NumberFieldElement_quadratic
63
import sage.rings.qqbar
64
AlgebraicNumber_base = sage.rings.qqbar.AlgebraicNumber_base
65
AlgebraicNumber = sage.rings.qqbar.AlgebraicNumber
66
AlgebraicReal = sage.rings.qqbar.AlgebraicReal
67
AA = sage.rings.qqbar.AA
68
QQbar = sage.rings.qqbar.QQbar
69
import sage.symbolic.ring
70
SR = sage.symbolic.ring.SR
71
from real_lazy import CLF, RLF
72
from complex_double import CDF
73
74
def is_ComplexField(x):
75
"""
76
Check if ``x`` is a :class:`complex field <ComplexField_class>`.
77
78
EXAMPLES::
79
80
sage: from sage.rings.complex_field import is_ComplexField as is_CF
81
sage: is_CF(ComplexField())
82
True
83
sage: is_CF(ComplexField(12))
84
True
85
sage: is_CF(CC)
86
True
87
"""
88
return isinstance(x, ComplexField_class)
89
90
cache = {}
91
def ComplexField(prec=53, names=None):
92
"""
93
Return the complex field with real and imaginary parts having prec
94
*bits* of precision.
95
96
EXAMPLES::
97
98
sage: ComplexField()
99
Complex Field with 53 bits of precision
100
sage: ComplexField(100)
101
Complex Field with 100 bits of precision
102
sage: ComplexField(100).base_ring()
103
Real Field with 100 bits of precision
104
sage: i = ComplexField(200).gen()
105
sage: i^2
106
-1.0000000000000000000000000000000000000000000000000000000000
107
"""
108
global cache
109
if cache.has_key(prec):
110
X = cache[prec]
111
C = X()
112
if not C is None:
113
return C
114
C = ComplexField_class(prec)
115
cache[prec] = weakref.ref(C)
116
return C
117
118
119
class ComplexField_class(field.Field):
120
"""
121
An approximation to the field of complex numbers using floating
122
point numbers with any specified precision. Answers derived from
123
calculations in this approximation may differ from what they would
124
be if those calculations were performed in the true field of
125
complex numbers. This is due to the rounding errors inherent to
126
finite precision calculations.
127
128
EXAMPLES::
129
130
sage: C = ComplexField(); C
131
Complex Field with 53 bits of precision
132
sage: Q = RationalField()
133
sage: C(1/3)
134
0.333333333333333
135
sage: C(1/3, 2)
136
0.333333333333333 + 2.00000000000000*I
137
sage: C(RR.pi())
138
3.14159265358979
139
sage: C(RR.log2(), RR.pi())
140
0.693147180559945 + 3.14159265358979*I
141
142
We can also coerce rational numbers and integers into C, but
143
coercing a polynomial will raise an exception::
144
145
sage: Q = RationalField()
146
sage: C(1/3)
147
0.333333333333333
148
sage: S = PolynomialRing(Q, 'x')
149
sage: C(S.gen())
150
Traceback (most recent call last):
151
...
152
TypeError: unable to coerce to a ComplexNumber: <type 'sage.rings.polynomial.polynomial_rational_flint.Polynomial_rational_flint'>
153
154
This illustrates precision::
155
156
sage: CC = ComplexField(10); CC(1/3, 2/3)
157
0.33 + 0.67*I
158
sage: CC
159
Complex Field with 10 bits of precision
160
sage: CC = ComplexField(100); CC
161
Complex Field with 100 bits of precision
162
sage: z = CC(1/3, 2/3); z
163
0.33333333333333333333333333333 + 0.66666666666666666666666666667*I
164
165
We can load and save complex numbers and the complex field::
166
167
sage: loads(z.dumps()) == z
168
True
169
sage: loads(CC.dumps()) == CC
170
True
171
sage: k = ComplexField(100)
172
sage: loads(dumps(k)) == k
173
True
174
175
This illustrates basic properties of a complex field::
176
177
sage: CC = ComplexField(200)
178
sage: CC.is_field()
179
True
180
sage: CC.characteristic()
181
0
182
sage: CC.precision()
183
200
184
sage: CC.variable_name()
185
'I'
186
sage: CC == ComplexField(200)
187
True
188
sage: CC == ComplexField(53)
189
False
190
sage: CC == 1.1
191
False
192
"""
193
def __init__(self, prec=53):
194
"""
195
Initialize ``self``.
196
197
TESTS::
198
199
sage: C = ComplexField(200)
200
sage: C.category()
201
Category of fields
202
sage: TestSuite(C).run()
203
"""
204
self._prec = int(prec)
205
from sage.categories.fields import Fields
206
ParentWithGens.__init__(self, self._real_field(), ('I',), False, category = Fields())
207
# self._populate_coercion_lists_()
208
self._populate_coercion_lists_(coerce_list=[complex_number.RRtoCC(self._real_field(), self)])
209
210
def __reduce__(self):
211
"""
212
For pickling.
213
214
EXAMPLES::
215
216
sage: loads(dumps(ComplexField())) == ComplexField()
217
True
218
"""
219
return ComplexField, (self._prec, )
220
221
def is_exact(self):
222
"""
223
Return whether or not this field is exact, which is always ``False``.
224
225
EXAMPLES::
226
227
sage: ComplexField().is_exact()
228
False
229
"""
230
return False
231
232
def prec(self):
233
"""
234
Return the precision of this complex field.
235
236
EXAMPLES::
237
238
sage: ComplexField().prec()
239
53
240
sage: ComplexField(15).prec()
241
15
242
"""
243
return self._prec
244
245
def _magma_init_(self, magma):
246
r"""
247
Return a string representation of ``self`` in the Magma language.
248
249
EXAMPLES::
250
251
sage: ComplexField()._magma_init_(magma) # optional - magma
252
'ComplexField(53 : Bits := true)'
253
sage: magma(ComplexField(200)) # optional - magma
254
Complex field of precision 60
255
sage: 10^60 < 2^200 < 10^61
256
True
257
sage: s = magma(ComplexField(200)).sage(); s # optional - magma
258
Complex Field with 200 bits of precision
259
sage: 2^199 < 10^60 < 2^200
260
True
261
sage: s is ComplexField(200) # optional - magma
262
True
263
"""
264
return "ComplexField(%s : Bits := true)" % self.prec()
265
266
precision = prec
267
268
def to_prec(self, prec):
269
"""
270
Returns the complex field to the specified precision.
271
272
EXAMPLES::
273
274
sage: CC.to_prec(10)
275
Complex Field with 10 bits of precision
276
sage: CC.to_prec(100)
277
Complex Field with 100 bits of precision
278
"""
279
return ComplexField(prec)
280
281
282
# very useful to cache this.
283
def _real_field(self):
284
"""
285
Return the underlying real field with the same precision.
286
287
EXAMPLES::
288
289
sage: RF = ComplexField(10)._real_field(); RF
290
Real Field with 10 bits of precision
291
sage: ComplexField(10)._real_field() is RF
292
True
293
"""
294
try:
295
return self.__real_field
296
except AttributeError:
297
self.__real_field = real_mpfr.RealField(self._prec)
298
return self.__real_field
299
300
def __cmp__(self, other):
301
"""
302
Compare ``self`` to ``other``.
303
304
If ``other`` is not a :class:`ComplexField_class', then this compares
305
by their types. Otherwise it compares by their precision.
306
307
EXAMPLES::
308
309
sage: cmp(ComplexField(), ComplexField())
310
0
311
sage: cmp(ComplexField(10), ComplexField(15))
312
-1
313
"""
314
if not isinstance(other, ComplexField_class):
315
return cmp(type(self), type(other))
316
return cmp(self._prec, other._prec)
317
318
def __call__(self, x=None, im=None):
319
"""
320
Create a complex number.
321
322
EXAMPLES::
323
324
sage: CC(2) # indirect doctest
325
2.00000000000000
326
sage: CC(CC.0)
327
1.00000000000000*I
328
sage: CC('1+I')
329
1.00000000000000 + 1.00000000000000*I
330
sage: CC(2,3)
331
2.00000000000000 + 3.00000000000000*I
332
sage: CC(QQ[I].gen())
333
1.00000000000000*I
334
sage: CC.gen() + QQ[I].gen()
335
Traceback (most recent call last):
336
...
337
TypeError: unsupported operand parent(s) for '+': 'Complex Field with 53 bits of precision' and 'Number Field in I with defining polynomial x^2 + 1'
338
339
In the absence of arguments we return zero::
340
341
sage: a = CC(); a
342
0.000000000000000
343
sage: a.parent()
344
Complex Field with 53 bits of precision
345
"""
346
if x is None:
347
return self.zero_element()
348
# we leave this here to handle the imaginary parameter
349
if im is not None:
350
x = x, im
351
return Parent.__call__(self, x)
352
353
def _element_constructor_(self, x):
354
"""
355
Construct a complex number.
356
357
EXAMPLES::
358
359
sage: CC((1,2)) # indirect doctest
360
1.00000000000000 + 2.00000000000000*I
361
362
Check that :trac:`14989` is fixed::
363
364
sage: QQi = NumberField(x^2+1, 'i', embedding=CC(0,1))
365
sage: i = QQi.order(QQi.gen()).gen(1)
366
sage: CC(i)
367
1.00000000000000*I
368
369
"""
370
if not isinstance(x, (real_mpfr.RealNumber, tuple)):
371
if isinstance(x, complex_double.ComplexDoubleElement):
372
return complex_number.ComplexNumber(self, x.real(), x.imag())
373
elif isinstance(x, str):
374
# TODO: this is probably not the best and most
375
# efficient way to do this. -- Martin Albrecht
376
return complex_number.ComplexNumber(self,
377
sage_eval(x.replace(' ',''), locals={"I":self.gen(),"i":self.gen()}))
378
379
late_import()
380
if isinstance(x, NumberFieldElement_quadratic):
381
if isinstance(x.parent(), NumberField_quadratic) and list(x.parent().polynomial()) == [1, 0, 1]:
382
(re, im) = list(x)
383
return complex_number.ComplexNumber(self, re, im)
384
385
try:
386
return self(x.sage())
387
except (AttributeError, TypeError):
388
pass
389
try:
390
return x._complex_mpfr_field_( self )
391
except AttributeError:
392
pass
393
return complex_number.ComplexNumber(self, x)
394
395
def _coerce_map_from_(self, S):
396
"""
397
The rings that canonically coerce to the MPFR complex field are:
398
399
- This MPFR complex field, or any other of higher precision
400
401
- Anything that canonically coerces to the mpfr real field
402
with this prec
403
404
EXAMPLES::
405
406
sage: ComplexField(200)(1) + RealField(90)(1) # indirect doctest
407
2.0000000000000000000000000
408
sage: parent(ComplexField(200)(1) + RealField(90)(1)) # indirect doctest
409
Complex Field with 90 bits of precision
410
sage: CC.0 + RLF(1/3) # indirect doctest
411
0.333333333333333 + 1.00000000000000*I
412
sage: ComplexField(20).has_coerce_map_from(CDF)
413
True
414
sage: ComplexField(200).has_coerce_map_from(CDF)
415
False
416
"""
417
RR = self._real_field()
418
if RR.has_coerce_map_from(S):
419
return complex_number.RRtoCC(RR, self) * RR.coerce_map_from(S)
420
if is_ComplexField(S) and S._prec >= self._prec:
421
return self._generic_convert_map(S)
422
late_import()
423
if S in [AA, QQbar, CLF, RLF] or (S == CDF and self._prec <= 53):
424
return self._generic_convert_map(S)
425
return self._coerce_map_via([CLF], S)
426
427
def _repr_(self):
428
"""
429
Return a string representation of ``self``.
430
431
EXAMPLES::
432
433
sage: ComplexField() # indirect doctest
434
Complex Field with 53 bits of precision
435
sage: ComplexField(15) # indirect doctest
436
Complex Field with 15 bits of precision
437
"""
438
return "Complex Field with %s bits of precision"%self._prec
439
440
def _latex_(self):
441
r"""
442
Return a latex representation of ``self``.
443
444
EXAMPLES::
445
446
sage: latex(ComplexField()) # indirect doctest
447
\Bold{C}
448
sage: latex(ComplexField(15)) # indirect doctest
449
\Bold{C}
450
"""
451
return "\\Bold{C}"
452
453
def _sage_input_(self, sib, coerce):
454
r"""
455
Produce an expression which will reproduce this value when evaluated.
456
457
EXAMPLES::
458
459
sage: sage_input(CC, verify=True)
460
# Verified
461
CC
462
sage: sage_input(ComplexField(25), verify=True)
463
# Verified
464
ComplexField(25)
465
sage: k = (CC, ComplexField(75))
466
sage: sage_input(k, verify=True)
467
# Verified
468
(CC, ComplexField(75))
469
sage: sage_input((k, k), verify=True)
470
# Verified
471
CC75 = ComplexField(75)
472
((CC, CC75), (CC, CC75))
473
sage: from sage.misc.sage_input import SageInputBuilder
474
sage: ComplexField(99)._sage_input_(SageInputBuilder(), False)
475
{call: {atomic:ComplexField}({atomic:99})}
476
"""
477
if self.prec() == 53:
478
return sib.name('CC')
479
480
v = sib.name('ComplexField')(sib.int(self.prec()))
481
482
name = 'CC%d' % (self.prec())
483
sib.cache(self, v, name)
484
return v
485
486
def characteristic(self):
487
r"""
488
Return the characteristic of `\CC`, which is 0.
489
490
EXAMPLES::
491
492
sage: ComplexField().characteristic()
493
0
494
"""
495
return integer.Integer(0)
496
497
def gen(self, n=0):
498
"""
499
Return the generator of the complex field.
500
501
EXAMPLES::
502
503
sage: ComplexField().gen(0)
504
1.00000000000000*I
505
"""
506
if n != 0:
507
raise IndexError, "n must be 0"
508
return complex_number.ComplexNumber(self, 0, 1)
509
510
def is_field(self, proof = True):
511
"""
512
Return ``True`` since the complex numbers are a field.
513
514
EXAMPLES::
515
516
sage: CC.is_field()
517
True
518
"""
519
return True
520
521
def is_finite(self):
522
"""
523
Return ``False`` since there are infinite number of complex numbers.
524
525
EXAMPLES::
526
527
sage: CC.is_finite()
528
False
529
"""
530
return False
531
532
def construction(self):
533
"""
534
Returns the functorial construction of ``self``, namely the algebraic
535
closure of the real field with the same precision.
536
537
EXAMPLES::
538
539
sage: c, S = CC.construction(); S
540
Real Field with 53 bits of precision
541
sage: CC == c(S)
542
True
543
"""
544
from sage.categories.pushout import AlgebraicClosureFunctor
545
return (AlgebraicClosureFunctor(), self._real_field())
546
547
548
def random_element(self, component_max=1, *args, **kwds):
549
r"""
550
Returns a uniformly distributed random number inside a square
551
centered on the origin (by default, the square `[-1,1] \times [-1,1]`).
552
553
Passes additional arguments and keywords to underlying real field.
554
555
EXAMPLES::
556
557
sage: [CC.random_element() for _ in range(5)]
558
[0.153636193785613 - 0.502987375247518*I,
559
0.609589964322241 - 0.948854594338216*I,
560
0.968393085385764 - 0.148483595843485*I,
561
-0.908976099636549 + 0.126219184235123*I,
562
0.461226845462901 - 0.0420335212948924*I]
563
sage: CC6 = ComplexField(6)
564
sage: [CC6.random_element(2^-20) for _ in range(5)]
565
[-5.4e-7 - 3.3e-7*I, 2.1e-7 + 8.0e-7*I, -4.8e-7 - 8.6e-7*I, -6.0e-8 + 2.7e-7*I, 6.0e-8 + 1.8e-7*I]
566
sage: [CC6.random_element(pi^20) for _ in range(5)]
567
[6.7e8 - 5.4e8*I, -9.4e8 + 5.0e9*I, 1.2e9 - 2.7e8*I, -2.3e9 - 4.0e9*I, 7.7e9 + 1.2e9*I]
568
569
Passes extra positional or keyword arguments through::
570
571
sage: [CC.random_element(distribution='1/n') for _ in range(5)]
572
[-0.900931453455899 - 0.932172283929307*I,
573
0.327862582226912 + 0.828104487111727*I,
574
0.246299162813240 + 0.588214960163442*I,
575
0.892970599589521 - 0.266744694790704*I,
576
0.878458776600692 - 0.905641181799996*I]
577
"""
578
size = self._real_field()(component_max)
579
re = self._real_field().random_element(-size, size, *args, **kwds)
580
im = self._real_field().random_element(-size, size, *args, **kwds)
581
return self(re, im)
582
583
def pi(self):
584
r"""
585
Returns `\pi` as a complex number.
586
587
EXAMPLES::
588
589
sage: ComplexField().pi()
590
3.14159265358979
591
sage: ComplexField(100).pi()
592
3.1415926535897932384626433833
593
"""
594
return self(self._real_field().pi())
595
596
def ngens(self):
597
r"""
598
The number of generators of this complex field as an `\RR`-algebra.
599
600
There is one generator, namely ``sqrt(-1)``.
601
602
EXAMPLES::
603
604
sage: ComplexField().ngens()
605
1
606
"""
607
return 1
608
609
def zeta(self, n=2):
610
"""
611
Return a primitive `n`-th root of unity.
612
613
INPUT:
614
615
- ``n`` - an integer (default: 2)
616
617
OUTPUT: a complex `n`-th root of unity.
618
619
EXAMPLES::
620
621
sage: C = ComplexField()
622
sage: C.zeta(2)
623
-1.00000000000000
624
sage: C.zeta(5)
625
0.309016994374947 + 0.951056516295154*I
626
"""
627
from integer import Integer
628
n = Integer(n)
629
if n == 1:
630
x = self(1)
631
elif n == 2:
632
x = self(-1)
633
elif n >= 3:
634
# Use De Moivre
635
# e^(2*pi*i/n) = cos(2pi/n) + i *sin(2pi/n)
636
RR = self._real_field()
637
pi = RR.pi()
638
z = 2*pi/n
639
x = complex_number.ComplexNumber(self, z.cos(), z.sin())
640
x._set_multiplicative_order( n )
641
return x
642
643
def scientific_notation(self, status=None):
644
"""
645
Set or return the scientific notation printing flag.
646
647
If this flag is ``True`` then complex numbers with this space as parent
648
print using scientific notation.
649
650
EXAMPLES::
651
652
sage: C = ComplexField()
653
sage: C((0.025, 2))
654
0.0250000000000000 + 2.00000000000000*I
655
sage: C.scientific_notation(True)
656
sage: C((0.025, 2))
657
2.50000000000000e-2 + 2.00000000000000e0*I
658
sage: C.scientific_notation(False)
659
sage: C((0.025, 2))
660
0.0250000000000000 + 2.00000000000000*I
661
"""
662
return self._real_field().scientific_notation(status)
663
664
def algebraic_closure(self):
665
"""
666
Return the algebraic closure of ``self`` (which is itself).
667
668
EXAMPLES::
669
670
sage: CC
671
Complex Field with 53 bits of precision
672
sage: CC.algebraic_closure()
673
Complex Field with 53 bits of precision
674
sage: CC = ComplexField(1000)
675
sage: CC.algebraic_closure() is CC
676
True
677
"""
678
return self
679
680
681