Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagelib
Path: blob/master/sage/rings/fraction_field_element.pyx
4069 views
1
"""
2
Fraction Field Elements
3
4
AUTHORS:
5
6
- William Stein (input from David Joyner, David Kohel, and Joe Wetherell)
7
8
- Sebastian Pancratz (2010-01-06): Rewrite of addition, multiplication and
9
derivative to use Henrici's algorithms [Ho72]
10
11
REFERENCES:
12
13
- [Ho72] E. Horowitz, "Algorithms for Rational Function Arithmetic
14
Operations", Annual ACM Symposium on Theory of Computing, Proceedings of
15
the Fourth Annual ACM Symposium on Theory of Computing, pp. 108--118, 1972
16
17
"""
18
19
#*****************************************************************************
20
#
21
# Sage: System for Algebra and Geometry Experimentation
22
#
23
# Copyright (C) 2005 William Stein <[email protected]>
24
#
25
# Distributed under the terms of the GNU General Public License (GPL)
26
#
27
# This code is distributed in the hope that it will be useful,
28
# but WITHOUT ANY WARRANTY; without even the implied warranty of
29
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
30
# General Public License for more details.
31
#
32
# The full text of the GPL is available at:
33
#
34
# http://www.gnu.org/licenses/
35
#*****************************************************************************
36
37
import operator
38
39
from sage.structure.element cimport FieldElement, ModuleElement, RingElement, \
40
Element
41
from sage.structure.element import parent
42
43
import integer_ring
44
from integer_ring import ZZ
45
from rational_field import QQ
46
47
import sage.misc.latex as latex
48
49
def is_FractionFieldElement(x):
50
"""
51
Returns whether or not x is of type FractionFieldElement.
52
53
EXAMPLES::
54
55
sage: from sage.rings.fraction_field_element import is_FractionFieldElement
56
sage: R.<x> = ZZ[]
57
sage: is_FractionFieldElement(x/2)
58
False
59
sage: is_FractionFieldElement(2/x)
60
True
61
sage: is_FractionFieldElement(1/3)
62
False
63
"""
64
return isinstance(x, FractionFieldElement)
65
66
cdef class FractionFieldElement(FieldElement):
67
"""
68
EXAMPLES::
69
70
sage: K, x = FractionField(PolynomialRing(QQ, 'x')).objgen()
71
sage: K
72
Fraction Field of Univariate Polynomial Ring in x over Rational Field
73
sage: loads(K.dumps()) == K
74
True
75
sage: f = (x^3 + x)/(17 - x^19); f
76
(x^3 + x)/(-x^19 + 17)
77
sage: loads(f.dumps()) == f
78
True
79
80
TESTS:
81
82
Test if #5451 is fixed::
83
84
sage: A = FiniteField(9,'theta')['t']
85
sage: K.<t> = FractionField(A)
86
sage: f= 2/(t^2+2*t); g =t^9/(t^18 + t^10 + t^2);f+g
87
(2*t^15 + 2*t^14 + 2*t^13 + 2*t^12 + 2*t^11 + 2*t^10 + 2*t^9 + t^7 + t^6 + t^5 + t^4 + t^3 + t^2 + t + 1)/(t^17 + t^9 + t)
88
89
Test if #8671 is fixed::
90
91
sage: P.<n> = QQ[]
92
sage: F = P.fraction_field()
93
sage: P.one_element()//F.one_element()
94
1
95
sage: F.one_element().quo_rem(F.one_element())
96
(1, 0)
97
"""
98
cdef object __numerator
99
cdef object __denominator
100
101
def __init__(self, parent, numerator, denominator=1,
102
coerce=True, reduce=True):
103
"""
104
EXAMPLES::
105
106
sage: from sage.rings.fraction_field_element import FractionFieldElement
107
sage: K.<x> = Frac(ZZ['x'])
108
sage: FractionFieldElement(K, x, 4)
109
x/4
110
sage: FractionFieldElement(K, x, x, reduce=False)
111
x/x
112
sage: f = FractionFieldElement(K, 'hi', 1, coerce=False, reduce=False)
113
sage: f.numerator()
114
'hi'
115
116
sage: x = var('x')
117
sage: K((x + 1)/(x^2 + x + 1))
118
(x + 1)/(x^2 + x + 1)
119
sage: K(355/113)
120
355/113
121
122
The next example failed before #4376::
123
124
sage: K(pari((x + 1)/(x^2 + x + 1)))
125
(x + 1)/(x^2 + x + 1)
126
127
"""
128
FieldElement.__init__(self, parent)
129
if coerce:
130
try:
131
self.__numerator = parent.ring()(numerator)
132
self.__denominator = parent.ring()(denominator)
133
except (TypeError, ValueError):
134
# workaround for symbolic ring
135
if denominator == 1 and hasattr(numerator, 'numerator'):
136
self.__numerator = parent.ring()(numerator.numerator())
137
self.__denominator = parent.ring()(numerator.denominator())
138
else:
139
raise
140
else:
141
self.__numerator = numerator
142
self.__denominator = denominator
143
if reduce and parent.is_exact():
144
try:
145
self.reduce()
146
except ArithmeticError:
147
pass
148
if self.__denominator.is_zero():
149
raise ZeroDivisionError, "fraction field element division by zero"
150
151
def _im_gens_(self, codomain, im_gens):
152
"""
153
EXAMPLES::
154
155
sage: F = ZZ['x,y'].fraction_field()
156
sage: x,y = F.gens()
157
sage: K = GF(7)['a,b'].fraction_field()
158
sage: a,b = K.gens()
159
160
::
161
162
sage: phi = F.hom([a+b, a*b], K)
163
sage: phi(x+y) # indirect doctest
164
a*b + a + b
165
166
::
167
168
sage: (x^2/y)._im_gens_(K, [a+b, a*b])
169
(a^2 + 2*a*b + b^2)/(a*b)
170
sage: (x^2/y)._im_gens_(K, [a, a*b])
171
a/b
172
"""
173
nnum = codomain.coerce(self.__numerator._im_gens_(codomain, im_gens))
174
nden = codomain.coerce(self.__denominator._im_gens_(codomain, im_gens))
175
return codomain.coerce(nnum/nden)
176
177
def reduce(self):
178
"""
179
Divides out the gcd of the numerator and denominator.
180
181
Automatically called for exact rings, but because it may be
182
numerically unstable for inexact rings it must be called manually
183
in that case.
184
185
EXAMPLES::
186
187
sage: R.<x> = RealField(10)[]
188
sage: f = (x^2+2*x+1)/(x+1); f
189
(x^2 + 2.0*x + 1.0)/(x + 1.0)
190
sage: f.reduce(); f
191
x + 1.0
192
"""
193
try:
194
g = self.__numerator.gcd(self.__denominator)
195
if not g.is_unit():
196
num, _ = self.__numerator.quo_rem(g)
197
den, _ = self.__denominator.quo_rem(g)
198
else:
199
num = self.__numerator
200
den = self.__denominator
201
if not den.is_one() and den.is_unit():
202
try:
203
num *= den.inverse_of_unit()
204
den = den.parent().one_element()
205
except:
206
pass
207
self.__numerator = num
208
self.__denominator = den
209
except AttributeError:
210
raise ArithmeticError, "unable to reduce because lack of gcd or quo_rem algorithm"
211
except TypeError:
212
raise ArithmeticError, "unable to reduce because gcd algorithm doesn't work on input"
213
except NotImplementedError:
214
raise ArithmeticError, "unable to reduce because gcd algorithm not implemented on input"
215
216
def __copy__(self):
217
"""
218
EXAMPLES::
219
220
sage: R.<x,y> = ZZ[]
221
sage: f = x/y+1; f
222
(x + y)/y
223
sage: copy(f)
224
(x + y)/y
225
"""
226
return self.__class__(self._parent, self.__numerator,
227
self.__denominator, coerce=False, reduce=False)
228
229
def numerator(self):
230
"""
231
EXAMPLES::
232
233
sage: R.<x,y> = ZZ[]
234
sage: f = x/y+1; f
235
(x + y)/y
236
sage: f.numerator()
237
x + y
238
"""
239
return self.__numerator
240
241
def denominator(self):
242
"""
243
EXAMPLES::
244
245
sage: R.<x,y> = ZZ[]
246
sage: f = x/y+1; f
247
(x + y)/y
248
sage: f.denominator()
249
y
250
"""
251
return self.__denominator
252
253
254
def is_square(self,root=False):
255
"""
256
Returns whether or not self is a perfect square. If the optional
257
argument root is True, then also returns a square root (or None,
258
if the fraction field element is not square).
259
260
INPUT:
261
262
263
- ``root`` - whether or not to also return a square
264
root (default: False)
265
266
267
OUTPUT:
268
269
270
- ``bool`` - whether or not a square
271
272
- ``object`` - (optional) an actual square root if
273
found, and None otherwise.
274
275
EXAMPLES::
276
277
sage: R.<t> = QQ[]
278
sage: (1/t).is_square()
279
False
280
sage: (1/t^6).is_square()
281
True
282
sage: ((1+t)^4/t^6).is_square()
283
True
284
sage: (4*(1+t)^4/t^6).is_square()
285
True
286
sage: (2*(1+t)^4/t^6).is_square()
287
False
288
sage: ((1+t)/t^6).is_square()
289
False
290
291
sage: (4*(1+t)^4/t^6).is_square(root=True)
292
(True, (2*t^2 + 4*t + 2)/t^3)
293
sage: (2*(1+t)^4/t^6).is_square(root=True)
294
(False, None)
295
296
sage: R.<x> = QQ[]
297
sage: a = 2*(x+1)^2 / (2*(x-1)^2); a
298
(2*x^2 + 4*x + 2)/(2*x^2 - 4*x + 2)
299
sage: a.numerator().is_square()
300
False
301
sage: a.is_square()
302
True
303
sage: (0/x).is_square()
304
True
305
"""
306
a = self.numerator()
307
b = self.denominator()
308
if not root:
309
return (a*b).is_square( root = False )
310
is_sqr, sq_rt = (a*b).is_square( root = True )
311
if is_sqr:
312
return True, self._parent( sq_rt/b )
313
return False, None
314
315
316
def __hash__(self):
317
"""
318
This function hashes in a special way to ensure that generators of
319
a ring `R` and generators of a fraction field of `R` have the same
320
hash. This enables them to be used as keys interchangeably in a
321
dictionary (since ``==`` will claim them equal). This is particularly
322
useful for methods like ``subs`` on ``ParentWithGens`` if you are
323
passing a dictionary of substitutions.
324
325
EXAMPLES::
326
327
sage: R.<x>=ZZ[]
328
sage: hash(R.0)==hash(FractionField(R).0)
329
True
330
sage: ((x+1)/(x^2+1)).subs({x:1})
331
1
332
sage: d={x:1}
333
sage: d[FractionField(R).0]
334
1
335
sage: R.<x>=QQ[] # this probably has a separate implementation from ZZ[]
336
sage: hash(R.0)==hash(FractionField(R).0)
337
True
338
sage: d={x:1}
339
sage: d[FractionField(R).0]
340
1
341
sage: R.<x,y,z>=ZZ[] # this probably has a separate implementation from ZZ[]
342
sage: hash(R.0)==hash(FractionField(R).0)
343
True
344
sage: d={x:1}
345
sage: d[FractionField(R).0]
346
1
347
sage: R.<x,y,z>=QQ[] # this probably has a separate implementation from ZZ[]
348
sage: hash(R.0)==hash(FractionField(R).0)
349
True
350
sage: ((x+1)/(x^2+1)).subs({x:1})
351
1
352
sage: d={x:1}
353
sage: d[FractionField(R).0]
354
1
355
sage: hash(R(1)/R(2))==hash(1/2)
356
True
357
"""
358
# This is same algorithm as used for members of QQ
359
#cdef long n, d
360
n = hash(self.__numerator)
361
d = hash(self.__denominator)
362
if d == 1:
363
return n
364
n = n ^ d
365
if n == -1:
366
return -2
367
return n
368
369
def __call__(self, *x, **kwds):
370
"""
371
Evaluate the fraction at the given arguments. This assumes that a
372
call function is defined for the numerator and denominator.
373
374
EXAMPLES::
375
376
sage: x = PolynomialRing(RationalField(),'x',3).gens()
377
sage: f = x[0] + x[1] - 2*x[1]*x[2]
378
sage: f
379
-2*x1*x2 + x0 + x1
380
sage: f(1,2,5)
381
-17
382
sage: h = f /(x[1] + x[2])
383
sage: h
384
(-2*x1*x2 + x0 + x1)/(x1 + x2)
385
sage: h(1,2,5)
386
-17/7
387
sage: h(x0=1)
388
(-2*x1*x2 + x1 + 1)/(x1 + x2)
389
"""
390
return self.__numerator(*x, **kwds) / self.__denominator(*x, **kwds)
391
392
def _is_atomic(self):
393
"""
394
EXAMPLES::
395
396
sage: K.<x> = Frac(ZZ['x'])
397
sage: x._is_atomic()
398
True
399
sage: f = 1/(x+1)
400
sage: f._is_atomic()
401
False
402
"""
403
return self.__numerator._is_atomic() and self.__denominator._is_atomic()
404
405
def _repr_(self):
406
"""
407
EXAMPLES::
408
409
sage: K.<x> = Frac(ZZ['x'])
410
sage: repr(x+1)
411
'x + 1'
412
sage: repr((x+1)/(x-1))
413
'(x + 1)/(x - 1)'
414
sage: repr(1/(x-1))
415
'1/(x - 1)'
416
sage: repr(1/x)
417
'1/x'
418
"""
419
if self.is_zero():
420
return "0"
421
s = "%s"%self.__numerator
422
if self.__denominator != 1:
423
denom_string = str( self.__denominator )
424
if self.__denominator._is_atomic() and not ('*' in denom_string or '/' in denom_string):
425
s = "%s/%s"%(self.__numerator._coeff_repr(no_space=False),denom_string)
426
else:
427
s = "%s/(%s)"%(self.__numerator._coeff_repr(no_space=False),denom_string)
428
return s
429
430
def _latex_(self):
431
r"""
432
Return a latex representation of this fraction field element.
433
434
EXAMPLES::
435
436
sage: R = PolynomialRing(QQ, 'x')
437
sage: F = R.fraction_field()
438
sage: x = F.gen()
439
sage: a = x^2 / 1
440
sage: latex(a)
441
x^{2}
442
sage: latex(x^2/(x^2+1))
443
\frac{x^{2}}{x^{2} + 1}
444
sage: a = 1/x
445
sage: latex(a)
446
\frac{1}{x}
447
448
TESTS::
449
450
sage: R = RR['x'] # Inexact, so no reduction.
451
sage: F = Frac(R)
452
sage: from sage.rings.fraction_field_element import FractionFieldElement
453
sage: z = FractionFieldElement(F, 0, R.gen(), coerce=False)
454
sage: z.numerator() == 0
455
True
456
sage: z.denominator() == R.gen()
457
True
458
sage: latex(z)
459
0
460
"""
461
if self.is_zero():
462
return "0"
463
if self.__denominator == 1:
464
return latex.latex(self.__numerator)
465
return "\\frac{%s}{%s}"%(latex.latex(self.__numerator),
466
latex.latex(self.__denominator))
467
468
def _magma_init_(self, magma):
469
"""
470
Return a string representation of ``self`` Magma can understand.
471
472
EXAMPLES::
473
474
sage: R.<x> = ZZ[]
475
sage: magma((x^2 + x + 1)/(x + 1)) # optional - magma
476
(x^2 + x + 1)/(x + 1)
477
478
::
479
480
sage: R.<x,y> = QQ[]
481
sage: magma((x+y)/x) # optional - magma
482
(x + y)/x
483
"""
484
pgens = magma(self._parent).gens()
485
486
s = self._repr_()
487
for i, j in zip(self._parent.variable_names(), pgens):
488
s = s.replace(i, j.name())
489
490
return s
491
492
cpdef ModuleElement _add_(self, ModuleElement right):
493
"""
494
Computes the sum of ``self`` and ``right``.
495
496
INPUT:
497
498
- ``right`` - ModuleElement to add to ``self``
499
500
OUTPUT:
501
502
- Sum of ``self`` and ``right``
503
504
EXAMPLES::
505
506
sage: K.<x,y> = Frac(ZZ['x,y'])
507
sage: x+y
508
x + y
509
sage: 1/x + 1/y
510
(x + y)/(x*y)
511
sage: 1/x + 1/(x*y)
512
(y + 1)/(x*y)
513
sage: Frac(CDF['x']).gen() + 3
514
x + 3.0
515
"""
516
rnum = self.__numerator
517
rden = self.__denominator
518
snum = (<FractionFieldElement> right).__numerator
519
sden = (<FractionFieldElement> right).__denominator
520
521
if (rnum.is_zero()):
522
return <FractionFieldElement> right
523
if (snum.is_zero()):
524
return self
525
526
if self._parent.is_exact():
527
try:
528
d = rden.gcd(sden)
529
if d.is_unit():
530
return self.__class__(self._parent, rnum*sden + rden*snum,
531
rden*sden, coerce=False, reduce=False)
532
else:
533
rden = rden // d
534
sden = sden // d
535
tnum = rnum * sden + rden * snum
536
if tnum.is_zero():
537
return self.__class__(self._parent, tnum,
538
self._parent.ring().one_element(), coerce=False,
539
reduce=False)
540
else:
541
tden = self.__denominator * sden
542
e = tnum.gcd(d)
543
if not e.is_unit():
544
tnum = tnum // e
545
tden = tden // e
546
if not tden.is_one() and tden.is_unit():
547
try:
548
tnum = tnum * tden.inverse_of_unit()
549
tden = self._parent.ring().one_element()
550
except AttributeError:
551
pass
552
except NotImplementedError:
553
pass
554
return self.__class__(self._parent, tnum, tden,
555
coerce=False, reduce=False)
556
except AttributeError:
557
pass
558
except NotImplementedError:
559
pass
560
except TypeError:
561
pass
562
563
rnum = self.__numerator
564
rden = self.__denominator
565
snum = (<FractionFieldElement> right).__numerator
566
sden = (<FractionFieldElement> right).__denominator
567
568
return self.__class__(self._parent, rnum*sden + rden*snum, rden*sden,
569
coerce=False, reduce=False)
570
571
cpdef ModuleElement _sub_(self, ModuleElement right):
572
"""
573
Computes the difference of ``self`` and ``right``.
574
575
INPUT:
576
577
- ``right`` - ModuleElement to subtract from ``self``
578
579
OUTPUT:
580
581
- Difference of ``self`` and ``right``
582
583
EXAMPLES::
584
585
sage: K.<t> = Frac(GF(7)['t'])
586
sage: t - 1/t
587
(t^2 + 6)/t
588
"""
589
return self._add_(-right)
590
591
cpdef RingElement _mul_(self, RingElement right):
592
"""
593
Computes the product of ``self`` and ``right``.
594
595
INPUT:
596
597
- ``right`` - RingElement to multiply with ``self``
598
599
OUTPUT:
600
601
- Product of ``self`` and ``right``
602
603
EXAMPLES::
604
605
sage: K.<t> = Frac(GF(7)['t'])
606
sage: a = t/(1+t)
607
sage: b = 3/t
608
sage: a*b
609
3/(t + 1)
610
"""
611
rnum = self.__numerator
612
rden = self.__denominator
613
snum = (<FractionFieldElement> right).__numerator
614
sden = (<FractionFieldElement> right).__denominator
615
616
if (rnum.is_zero() or snum.is_zero()):
617
return self._parent.zero_element()
618
619
if self._parent.is_exact():
620
try:
621
d1 = rnum.gcd(sden)
622
d2 = snum.gcd(rden)
623
if not d1.is_unit():
624
rnum = rnum // d1
625
sden = sden // d1
626
if not d2.is_unit():
627
rden = rden // d2
628
snum = snum // d2
629
tnum = rnum * snum
630
tden = rden * sden
631
if not tden.is_one() and tden.is_unit():
632
try:
633
tnum = tnum * tden.inverse_of_unit()
634
tden = self._parent.ring().one_element()
635
except AttributeError:
636
pass
637
except NotImplementedError:
638
pass
639
return self.__class__(self._parent, tnum, tden,
640
coerce=False, reduce=False)
641
except AttributeError:
642
pass
643
except NotImplementedError:
644
pass
645
except TypeError:
646
pass
647
648
rnum = self.__numerator
649
rden = self.__denominator
650
snum = (<FractionFieldElement> right).__numerator
651
sden = (<FractionFieldElement> right).__denominator
652
653
return self.__class__(self._parent, rnum * snum, rden * sden,
654
coerce=False, reduce=False)
655
656
cpdef RingElement _div_(self, RingElement right):
657
"""
658
Computes the quotient of ``self`` and ``right``.
659
660
INPUT:
661
662
- ``right`` - RingElement that is the divisor
663
664
OUTPUT:
665
666
- Quotient of ``self`` and ``right``
667
668
EXAMPLES::
669
670
sage: K.<x,y,z> = Frac(ZZ['x,y,z'])
671
sage: a = (x+1)*(x+y)/(z-3)
672
sage: b = (x+y)/(z-1)
673
sage: a/b
674
(x*z - x + z - 1)/(z - 3)
675
"""
676
snum = (<FractionFieldElement> right).__numerator
677
sden = (<FractionFieldElement> right).__denominator
678
679
if snum.is_zero():
680
raise ZeroDivisionError, "fraction field element division by zero"
681
682
rightinv = self.__class__(self._parent, sden, snum,
683
coerce=True, reduce=False)
684
685
return self._mul_(rightinv)
686
687
def __int__(self):
688
"""
689
EXAMPLES::
690
691
sage: K = Frac(ZZ['x'])
692
sage: int(K(-3))
693
-3
694
sage: K.<x> = Frac(RR['x'])
695
sage: x/x
696
x/x
697
sage: int(x/x)
698
1
699
sage: int(K(.5))
700
0
701
"""
702
if self.__denominator != 1:
703
self.reduce()
704
if self.__denominator == 1:
705
return int(self.__numerator)
706
else:
707
raise TypeError, "denominator must equal 1"
708
709
def _integer_(self, Z=ZZ):
710
"""
711
EXAMPLES::
712
713
sage: K = Frac(ZZ['x'])
714
sage: K(5)._integer_()
715
5
716
sage: K.<x> = Frac(RR['x'])
717
sage: ZZ(2*x/x)
718
2
719
"""
720
if self.__denominator != 1:
721
self.reduce()
722
if self.__denominator == 1:
723
return Z(self.__numerator)
724
raise TypeError, "no way to coerce to an integer."
725
726
def _rational_(self, Q=QQ):
727
"""
728
EXAMPLES::
729
730
sage: K.<x> = Frac(QQ['x'])
731
sage: K(1/2)._rational_()
732
1/2
733
sage: K(1/2 + x/x)._rational_()
734
3/2
735
"""
736
return Q(self.__numerator) / Q(self.__denominator)
737
738
def __long__(self):
739
"""
740
EXAMPLES::
741
742
sage: K.<x> = Frac(QQ['x'])
743
sage: long(K(3))
744
3L
745
sage: long(K(3/5))
746
0L
747
"""
748
return long(int(self))
749
750
def __pow__(self, right, dummy):
751
r"""
752
Returns self raised to the `right^{th}` power.
753
754
Note that we need to check whether or not right is negative so we
755
don't set __numerator or __denominator to an element of the
756
fraction field instead of the underlying ring.
757
758
EXAMPLES::
759
760
sage: R = QQ['x','y']
761
sage: FR = R.fraction_field()
762
sage: x,y = FR.gens()
763
sage: a = x^2; a
764
x^2
765
sage: type(a.numerator())
766
<type 'sage.rings.polynomial.multi_polynomial_libsingular.MPolynomial_libsingular'>
767
sage: type(a.denominator())
768
<type 'sage.rings.polynomial.multi_polynomial_libsingular.MPolynomial_libsingular'>
769
sage: a = x^(-2); a
770
1/x^2
771
sage: type(a.numerator())
772
<type 'sage.rings.polynomial.multi_polynomial_libsingular.MPolynomial_libsingular'>
773
sage: type(a.denominator())
774
<type 'sage.rings.polynomial.multi_polynomial_libsingular.MPolynomial_libsingular'>
775
sage: x^0
776
1
777
sage: ((x+y)/(x-y))^2
778
(x^2 + 2*x*y + y^2)/(x^2 - 2*x*y + y^2)
779
sage: ((x+y)/(x-y))^-2
780
(x^2 - 2*x*y + y^2)/(x^2 + 2*x*y + y^2)
781
sage: ((x+y)/(x-y))^0
782
1
783
"""
784
snum = (<FractionFieldElement> self).__numerator
785
sden = (<FractionFieldElement> self).__denominator
786
if right == 0:
787
R = self.parent().ring()
788
return self.__class__(self.parent(),
789
R.one_element(), R.one_element(),
790
coerce=False, reduce=False)
791
elif right > 0:
792
return self.__class__(self.parent(),
793
snum**right, sden**right,
794
coerce=False, reduce=False)
795
else:
796
right = -right
797
return self.__class__(self.parent(),
798
sden**right, snum**right,
799
coerce=False, reduce=False)
800
801
def __neg__(self):
802
"""
803
EXAMPLES::
804
805
sage: K.<t> = Frac(GF(5)['t'])
806
sage: f = (t^2+t)/(t+2); f
807
(t^2 + t)/(t + 2)
808
sage: -f
809
(4*t^2 + 4*t)/(t + 2)
810
"""
811
return self.__class__(self._parent,
812
-self.__numerator, self.__denominator,
813
coerce=False, reduce=False)
814
815
def __abs__(self):
816
"""
817
EXAMPLES::
818
819
sage: from sage.rings.fraction_field_element import FractionFieldElement
820
sage: abs(FractionFieldElement(QQ, -2, 3, coerce=False))
821
2/3
822
"""
823
return abs(self.__numerator)/abs(self.__denominator)
824
825
def __invert__(self):
826
"""
827
EXAMPLES::
828
829
sage: K.<t> = Frac(GF(7)['t'])
830
sage: f = (t^2+5)/(t-1)
831
sage: ~f
832
(t + 6)/(t^2 + 5)
833
"""
834
if self.is_zero():
835
raise ZeroDivisionError, "Cannot invert 0"
836
return self.__class__(self._parent,
837
self.__denominator, self.__numerator, coerce=False, reduce=False)
838
839
def __float__(self):
840
"""
841
EXAMPLES::
842
843
sage: K.<x,y> = Frac(ZZ['x,y'])
844
sage: float(x/x + y/y)
845
2.0
846
"""
847
return float(self.__numerator) / float(self.__denominator)
848
849
def __richcmp__(left, right, int op):
850
"""
851
EXAMPLES::
852
853
sage: K.<x,y> = Frac(ZZ['x,y'])
854
sage: x > y
855
True
856
sage: 1 > y
857
False
858
"""
859
return (<Element>left)._richcmp(right, op)
860
861
cdef int _cmp_c_impl(self, Element other) except -2:
862
"""
863
EXAMPLES::
864
865
sage: K.<t> = Frac(GF(7)['t'])
866
sage: t/t == 1
867
True
868
sage: t+1/t == (t^2+1)/t
869
True
870
sage: t == t/5
871
False
872
"""
873
return cmp(self.__numerator * \
874
(<FractionFieldElement>other).__denominator,
875
self.__denominator*(<FractionFieldElement>other).__numerator)
876
877
def valuation(self, v=None):
878
"""
879
Return the valuation of self, assuming that the numerator and
880
denominator have valuation functions defined on them.
881
882
EXAMPLES::
883
884
sage: x = PolynomialRing(RationalField(),'x').gen()
885
sage: f = (x^3 + x)/(x^2 - 2*x^3)
886
sage: f
887
(x^2 + 1)/(-2*x^2 + x)
888
sage: f.valuation()
889
-1
890
sage: f.valuation(x^2+1)
891
1
892
"""
893
return self.__numerator.valuation(v) - self.__denominator.valuation(v)
894
895
def __nonzero__(self):
896
"""
897
Returns True if this element is nonzero.
898
899
EXAMPLES::
900
901
sage: F = ZZ['x,y'].fraction_field()
902
sage: x,y = F.gens()
903
sage: t = F(0)/x
904
sage: t.__nonzero__()
905
False
906
907
::
908
909
sage: (1/x).__nonzero__()
910
True
911
"""
912
return not self.__numerator.is_zero()
913
914
def is_zero(self):
915
"""
916
Returns True if this element is equal to zero.
917
918
EXAMPLES::
919
920
sage: F = ZZ['x,y'].fraction_field()
921
sage: x,y = F.gens()
922
sage: t = F(0)/x
923
sage: t.is_zero()
924
True
925
sage: u = 1/x - 1/x
926
sage: u.is_zero()
927
True
928
sage: u.parent() is F
929
True
930
"""
931
return self.__numerator.is_zero()
932
933
def is_one(self):
934
"""
935
Returns True if this element is equal to one.
936
937
EXAMPLES::
938
939
sage: F = ZZ['x,y'].fraction_field()
940
sage: x,y = F.gens()
941
sage: (x/x).is_one()
942
True
943
sage: (x/y).is_one()
944
False
945
"""
946
return self.__numerator == self.__denominator
947
948
def _symbolic_(self, ring):
949
return ring(self.__numerator)/ring(self.__denominator)
950
951
def __reduce__(self):
952
"""
953
EXAMPLES::
954
955
sage: F = ZZ['x,y'].fraction_field()
956
sage: f = F.random_element()
957
sage: loads(f.dumps()) == f
958
True
959
"""
960
return (make_element,
961
(self._parent, self.__numerator, self.__denominator))
962
963
964
class FractionFieldElement_1poly_field(FractionFieldElement):
965
"""
966
A fraction field element where the parent is the fraction field of a univariate polynomial ring.
967
968
Many of the functions here are included for coherence with number fields.
969
"""
970
def is_integral(self):
971
"""
972
Returns whether this element is actually a polynomial.
973
974
EXAMPLES::
975
976
sage: R.<t> = GF(5)[]
977
sage: K = R.fraction_field
978
"""
979
if self.__denominator != 1:
980
self.reduce()
981
return self.__denominator == 1
982
983
def support(self):
984
"""
985
Returns a sorted list of primes dividing either the numerator or denominator of this element.
986
987
EXAMPLES::
988
989
sage: R.<t> = QQ[]
990
sage: h = (t^14 + 2*t^12 - 4*t^11 - 8*t^9 + 6*t^8 + 12*t^6 - 4*t^5 - 8*t^3 + t^2 + 2)/(t^6 + 6*t^5 + 9*t^4 - 2*t^2 - 12*t - 18)
991
sage: h.support()
992
[t - 1, t + 3, t^2 + 2, t^2 + t + 1, t^4 - 2]
993
"""
994
L = [fac[0] for fac in self.numerator().factor()] + [fac[0] for fac in self.denominator().factor()]
995
L.sort()
996
return L
997
998
999
1000
def make_element(parent, numerator, denominator):
1001
"""
1002
Used for unpickling FractionFieldElement objects (and subclasses).
1003
1004
EXAMPLES::
1005
1006
sage: from sage.rings.fraction_field_element import make_element
1007
sage: R = ZZ['x,y']
1008
sage: x,y = R.gens()
1009
sage: F = R.fraction_field()
1010
sage: make_element(F, 1+x, 1+y)
1011
(x + 1)/(y + 1)
1012
"""
1013
1014
return parent._element_class(parent, numerator, denominator)
1015
1016
def make_element_old(parent, cdict):
1017
"""
1018
Used for unpickling old FractionFieldElement pickles.
1019
1020
EXAMPLES::
1021
1022
sage: from sage.rings.fraction_field_element import make_element_old
1023
sage: R.<x,y> = ZZ[]
1024
sage: F = R.fraction_field()
1025
sage: make_element_old(F, {'_FractionFieldElement__numerator':x+y,'_FractionFieldElement__denominator':x-y})
1026
(x + y)/(x - y)
1027
"""
1028
return FractionFieldElement(parent,
1029
cdict['_FractionFieldElement__numerator'],
1030
cdict['_FractionFieldElement__denominator'],
1031
coerce=False, reduce=False)
1032
1033
1034