Path: blob/master/sage/rings/fraction_field_element.pyx
4069 views
"""1Fraction Field Elements23AUTHORS:45- William Stein (input from David Joyner, David Kohel, and Joe Wetherell)67- Sebastian Pancratz (2010-01-06): Rewrite of addition, multiplication and8derivative to use Henrici's algorithms [Ho72]910REFERENCES:1112- [Ho72] E. Horowitz, "Algorithms for Rational Function Arithmetic13Operations", Annual ACM Symposium on Theory of Computing, Proceedings of14the Fourth Annual ACM Symposium on Theory of Computing, pp. 108--118, 19721516"""1718#*****************************************************************************19#20# Sage: System for Algebra and Geometry Experimentation21#22# Copyright (C) 2005 William Stein <[email protected]>23#24# Distributed under the terms of the GNU General Public License (GPL)25#26# This code is distributed in the hope that it will be useful,27# but WITHOUT ANY WARRANTY; without even the implied warranty of28# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU29# General Public License for more details.30#31# The full text of the GPL is available at:32#33# http://www.gnu.org/licenses/34#*****************************************************************************3536import operator3738from sage.structure.element cimport FieldElement, ModuleElement, RingElement, \39Element40from sage.structure.element import parent4142import integer_ring43from integer_ring import ZZ44from rational_field import QQ4546import sage.misc.latex as latex4748def is_FractionFieldElement(x):49"""50Returns whether or not x is of type FractionFieldElement.5152EXAMPLES::5354sage: from sage.rings.fraction_field_element import is_FractionFieldElement55sage: R.<x> = ZZ[]56sage: is_FractionFieldElement(x/2)57False58sage: is_FractionFieldElement(2/x)59True60sage: is_FractionFieldElement(1/3)61False62"""63return isinstance(x, FractionFieldElement)6465cdef class FractionFieldElement(FieldElement):66"""67EXAMPLES::6869sage: K, x = FractionField(PolynomialRing(QQ, 'x')).objgen()70sage: K71Fraction Field of Univariate Polynomial Ring in x over Rational Field72sage: loads(K.dumps()) == K73True74sage: f = (x^3 + x)/(17 - x^19); f75(x^3 + x)/(-x^19 + 17)76sage: loads(f.dumps()) == f77True7879TESTS:8081Test if #5451 is fixed::8283sage: A = FiniteField(9,'theta')['t']84sage: K.<t> = FractionField(A)85sage: f= 2/(t^2+2*t); g =t^9/(t^18 + t^10 + t^2);f+g86(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)8788Test if #8671 is fixed::8990sage: P.<n> = QQ[]91sage: F = P.fraction_field()92sage: P.one_element()//F.one_element()93194sage: F.one_element().quo_rem(F.one_element())95(1, 0)96"""97cdef object __numerator98cdef object __denominator99100def __init__(self, parent, numerator, denominator=1,101coerce=True, reduce=True):102"""103EXAMPLES::104105sage: from sage.rings.fraction_field_element import FractionFieldElement106sage: K.<x> = Frac(ZZ['x'])107sage: FractionFieldElement(K, x, 4)108x/4109sage: FractionFieldElement(K, x, x, reduce=False)110x/x111sage: f = FractionFieldElement(K, 'hi', 1, coerce=False, reduce=False)112sage: f.numerator()113'hi'114115sage: x = var('x')116sage: K((x + 1)/(x^2 + x + 1))117(x + 1)/(x^2 + x + 1)118sage: K(355/113)119355/113120121The next example failed before #4376::122123sage: K(pari((x + 1)/(x^2 + x + 1)))124(x + 1)/(x^2 + x + 1)125126"""127FieldElement.__init__(self, parent)128if coerce:129try:130self.__numerator = parent.ring()(numerator)131self.__denominator = parent.ring()(denominator)132except (TypeError, ValueError):133# workaround for symbolic ring134if denominator == 1 and hasattr(numerator, 'numerator'):135self.__numerator = parent.ring()(numerator.numerator())136self.__denominator = parent.ring()(numerator.denominator())137else:138raise139else:140self.__numerator = numerator141self.__denominator = denominator142if reduce and parent.is_exact():143try:144self.reduce()145except ArithmeticError:146pass147if self.__denominator.is_zero():148raise ZeroDivisionError, "fraction field element division by zero"149150def _im_gens_(self, codomain, im_gens):151"""152EXAMPLES::153154sage: F = ZZ['x,y'].fraction_field()155sage: x,y = F.gens()156sage: K = GF(7)['a,b'].fraction_field()157sage: a,b = K.gens()158159::160161sage: phi = F.hom([a+b, a*b], K)162sage: phi(x+y) # indirect doctest163a*b + a + b164165::166167sage: (x^2/y)._im_gens_(K, [a+b, a*b])168(a^2 + 2*a*b + b^2)/(a*b)169sage: (x^2/y)._im_gens_(K, [a, a*b])170a/b171"""172nnum = codomain.coerce(self.__numerator._im_gens_(codomain, im_gens))173nden = codomain.coerce(self.__denominator._im_gens_(codomain, im_gens))174return codomain.coerce(nnum/nden)175176def reduce(self):177"""178Divides out the gcd of the numerator and denominator.179180Automatically called for exact rings, but because it may be181numerically unstable for inexact rings it must be called manually182in that case.183184EXAMPLES::185186sage: R.<x> = RealField(10)[]187sage: f = (x^2+2*x+1)/(x+1); f188(x^2 + 2.0*x + 1.0)/(x + 1.0)189sage: f.reduce(); f190x + 1.0191"""192try:193g = self.__numerator.gcd(self.__denominator)194if not g.is_unit():195num, _ = self.__numerator.quo_rem(g)196den, _ = self.__denominator.quo_rem(g)197else:198num = self.__numerator199den = self.__denominator200if not den.is_one() and den.is_unit():201try:202num *= den.inverse_of_unit()203den = den.parent().one_element()204except:205pass206self.__numerator = num207self.__denominator = den208except AttributeError:209raise ArithmeticError, "unable to reduce because lack of gcd or quo_rem algorithm"210except TypeError:211raise ArithmeticError, "unable to reduce because gcd algorithm doesn't work on input"212except NotImplementedError:213raise ArithmeticError, "unable to reduce because gcd algorithm not implemented on input"214215def __copy__(self):216"""217EXAMPLES::218219sage: R.<x,y> = ZZ[]220sage: f = x/y+1; f221(x + y)/y222sage: copy(f)223(x + y)/y224"""225return self.__class__(self._parent, self.__numerator,226self.__denominator, coerce=False, reduce=False)227228def numerator(self):229"""230EXAMPLES::231232sage: R.<x,y> = ZZ[]233sage: f = x/y+1; f234(x + y)/y235sage: f.numerator()236x + y237"""238return self.__numerator239240def denominator(self):241"""242EXAMPLES::243244sage: R.<x,y> = ZZ[]245sage: f = x/y+1; f246(x + y)/y247sage: f.denominator()248y249"""250return self.__denominator251252253def is_square(self,root=False):254"""255Returns whether or not self is a perfect square. If the optional256argument root is True, then also returns a square root (or None,257if the fraction field element is not square).258259INPUT:260261262- ``root`` - whether or not to also return a square263root (default: False)264265266OUTPUT:267268269- ``bool`` - whether or not a square270271- ``object`` - (optional) an actual square root if272found, and None otherwise.273274EXAMPLES::275276sage: R.<t> = QQ[]277sage: (1/t).is_square()278False279sage: (1/t^6).is_square()280True281sage: ((1+t)^4/t^6).is_square()282True283sage: (4*(1+t)^4/t^6).is_square()284True285sage: (2*(1+t)^4/t^6).is_square()286False287sage: ((1+t)/t^6).is_square()288False289290sage: (4*(1+t)^4/t^6).is_square(root=True)291(True, (2*t^2 + 4*t + 2)/t^3)292sage: (2*(1+t)^4/t^6).is_square(root=True)293(False, None)294295sage: R.<x> = QQ[]296sage: a = 2*(x+1)^2 / (2*(x-1)^2); a297(2*x^2 + 4*x + 2)/(2*x^2 - 4*x + 2)298sage: a.numerator().is_square()299False300sage: a.is_square()301True302sage: (0/x).is_square()303True304"""305a = self.numerator()306b = self.denominator()307if not root:308return (a*b).is_square( root = False )309is_sqr, sq_rt = (a*b).is_square( root = True )310if is_sqr:311return True, self._parent( sq_rt/b )312return False, None313314315def __hash__(self):316"""317This function hashes in a special way to ensure that generators of318a ring `R` and generators of a fraction field of `R` have the same319hash. This enables them to be used as keys interchangeably in a320dictionary (since ``==`` will claim them equal). This is particularly321useful for methods like ``subs`` on ``ParentWithGens`` if you are322passing a dictionary of substitutions.323324EXAMPLES::325326sage: R.<x>=ZZ[]327sage: hash(R.0)==hash(FractionField(R).0)328True329sage: ((x+1)/(x^2+1)).subs({x:1})3301331sage: d={x:1}332sage: d[FractionField(R).0]3331334sage: R.<x>=QQ[] # this probably has a separate implementation from ZZ[]335sage: hash(R.0)==hash(FractionField(R).0)336True337sage: d={x:1}338sage: d[FractionField(R).0]3391340sage: R.<x,y,z>=ZZ[] # this probably has a separate implementation from ZZ[]341sage: hash(R.0)==hash(FractionField(R).0)342True343sage: d={x:1}344sage: d[FractionField(R).0]3451346sage: R.<x,y,z>=QQ[] # this probably has a separate implementation from ZZ[]347sage: hash(R.0)==hash(FractionField(R).0)348True349sage: ((x+1)/(x^2+1)).subs({x:1})3501351sage: d={x:1}352sage: d[FractionField(R).0]3531354sage: hash(R(1)/R(2))==hash(1/2)355True356"""357# This is same algorithm as used for members of QQ358#cdef long n, d359n = hash(self.__numerator)360d = hash(self.__denominator)361if d == 1:362return n363n = n ^ d364if n == -1:365return -2366return n367368def __call__(self, *x, **kwds):369"""370Evaluate the fraction at the given arguments. This assumes that a371call function is defined for the numerator and denominator.372373EXAMPLES::374375sage: x = PolynomialRing(RationalField(),'x',3).gens()376sage: f = x[0] + x[1] - 2*x[1]*x[2]377sage: f378-2*x1*x2 + x0 + x1379sage: f(1,2,5)380-17381sage: h = f /(x[1] + x[2])382sage: h383(-2*x1*x2 + x0 + x1)/(x1 + x2)384sage: h(1,2,5)385-17/7386sage: h(x0=1)387(-2*x1*x2 + x1 + 1)/(x1 + x2)388"""389return self.__numerator(*x, **kwds) / self.__denominator(*x, **kwds)390391def _is_atomic(self):392"""393EXAMPLES::394395sage: K.<x> = Frac(ZZ['x'])396sage: x._is_atomic()397True398sage: f = 1/(x+1)399sage: f._is_atomic()400False401"""402return self.__numerator._is_atomic() and self.__denominator._is_atomic()403404def _repr_(self):405"""406EXAMPLES::407408sage: K.<x> = Frac(ZZ['x'])409sage: repr(x+1)410'x + 1'411sage: repr((x+1)/(x-1))412'(x + 1)/(x - 1)'413sage: repr(1/(x-1))414'1/(x - 1)'415sage: repr(1/x)416'1/x'417"""418if self.is_zero():419return "0"420s = "%s"%self.__numerator421if self.__denominator != 1:422denom_string = str( self.__denominator )423if self.__denominator._is_atomic() and not ('*' in denom_string or '/' in denom_string):424s = "%s/%s"%(self.__numerator._coeff_repr(no_space=False),denom_string)425else:426s = "%s/(%s)"%(self.__numerator._coeff_repr(no_space=False),denom_string)427return s428429def _latex_(self):430r"""431Return a latex representation of this fraction field element.432433EXAMPLES::434435sage: R = PolynomialRing(QQ, 'x')436sage: F = R.fraction_field()437sage: x = F.gen()438sage: a = x^2 / 1439sage: latex(a)440x^{2}441sage: latex(x^2/(x^2+1))442\frac{x^{2}}{x^{2} + 1}443sage: a = 1/x444sage: latex(a)445\frac{1}{x}446447TESTS::448449sage: R = RR['x'] # Inexact, so no reduction.450sage: F = Frac(R)451sage: from sage.rings.fraction_field_element import FractionFieldElement452sage: z = FractionFieldElement(F, 0, R.gen(), coerce=False)453sage: z.numerator() == 0454True455sage: z.denominator() == R.gen()456True457sage: latex(z)4580459"""460if self.is_zero():461return "0"462if self.__denominator == 1:463return latex.latex(self.__numerator)464return "\\frac{%s}{%s}"%(latex.latex(self.__numerator),465latex.latex(self.__denominator))466467def _magma_init_(self, magma):468"""469Return a string representation of ``self`` Magma can understand.470471EXAMPLES::472473sage: R.<x> = ZZ[]474sage: magma((x^2 + x + 1)/(x + 1)) # optional - magma475(x^2 + x + 1)/(x + 1)476477::478479sage: R.<x,y> = QQ[]480sage: magma((x+y)/x) # optional - magma481(x + y)/x482"""483pgens = magma(self._parent).gens()484485s = self._repr_()486for i, j in zip(self._parent.variable_names(), pgens):487s = s.replace(i, j.name())488489return s490491cpdef ModuleElement _add_(self, ModuleElement right):492"""493Computes the sum of ``self`` and ``right``.494495INPUT:496497- ``right`` - ModuleElement to add to ``self``498499OUTPUT:500501- Sum of ``self`` and ``right``502503EXAMPLES::504505sage: K.<x,y> = Frac(ZZ['x,y'])506sage: x+y507x + y508sage: 1/x + 1/y509(x + y)/(x*y)510sage: 1/x + 1/(x*y)511(y + 1)/(x*y)512sage: Frac(CDF['x']).gen() + 3513x + 3.0514"""515rnum = self.__numerator516rden = self.__denominator517snum = (<FractionFieldElement> right).__numerator518sden = (<FractionFieldElement> right).__denominator519520if (rnum.is_zero()):521return <FractionFieldElement> right522if (snum.is_zero()):523return self524525if self._parent.is_exact():526try:527d = rden.gcd(sden)528if d.is_unit():529return self.__class__(self._parent, rnum*sden + rden*snum,530rden*sden, coerce=False, reduce=False)531else:532rden = rden // d533sden = sden // d534tnum = rnum * sden + rden * snum535if tnum.is_zero():536return self.__class__(self._parent, tnum,537self._parent.ring().one_element(), coerce=False,538reduce=False)539else:540tden = self.__denominator * sden541e = tnum.gcd(d)542if not e.is_unit():543tnum = tnum // e544tden = tden // e545if not tden.is_one() and tden.is_unit():546try:547tnum = tnum * tden.inverse_of_unit()548tden = self._parent.ring().one_element()549except AttributeError:550pass551except NotImplementedError:552pass553return self.__class__(self._parent, tnum, tden,554coerce=False, reduce=False)555except AttributeError:556pass557except NotImplementedError:558pass559except TypeError:560pass561562rnum = self.__numerator563rden = self.__denominator564snum = (<FractionFieldElement> right).__numerator565sden = (<FractionFieldElement> right).__denominator566567return self.__class__(self._parent, rnum*sden + rden*snum, rden*sden,568coerce=False, reduce=False)569570cpdef ModuleElement _sub_(self, ModuleElement right):571"""572Computes the difference of ``self`` and ``right``.573574INPUT:575576- ``right`` - ModuleElement to subtract from ``self``577578OUTPUT:579580- Difference of ``self`` and ``right``581582EXAMPLES::583584sage: K.<t> = Frac(GF(7)['t'])585sage: t - 1/t586(t^2 + 6)/t587"""588return self._add_(-right)589590cpdef RingElement _mul_(self, RingElement right):591"""592Computes the product of ``self`` and ``right``.593594INPUT:595596- ``right`` - RingElement to multiply with ``self``597598OUTPUT:599600- Product of ``self`` and ``right``601602EXAMPLES::603604sage: K.<t> = Frac(GF(7)['t'])605sage: a = t/(1+t)606sage: b = 3/t607sage: a*b6083/(t + 1)609"""610rnum = self.__numerator611rden = self.__denominator612snum = (<FractionFieldElement> right).__numerator613sden = (<FractionFieldElement> right).__denominator614615if (rnum.is_zero() or snum.is_zero()):616return self._parent.zero_element()617618if self._parent.is_exact():619try:620d1 = rnum.gcd(sden)621d2 = snum.gcd(rden)622if not d1.is_unit():623rnum = rnum // d1624sden = sden // d1625if not d2.is_unit():626rden = rden // d2627snum = snum // d2628tnum = rnum * snum629tden = rden * sden630if not tden.is_one() and tden.is_unit():631try:632tnum = tnum * tden.inverse_of_unit()633tden = self._parent.ring().one_element()634except AttributeError:635pass636except NotImplementedError:637pass638return self.__class__(self._parent, tnum, tden,639coerce=False, reduce=False)640except AttributeError:641pass642except NotImplementedError:643pass644except TypeError:645pass646647rnum = self.__numerator648rden = self.__denominator649snum = (<FractionFieldElement> right).__numerator650sden = (<FractionFieldElement> right).__denominator651652return self.__class__(self._parent, rnum * snum, rden * sden,653coerce=False, reduce=False)654655cpdef RingElement _div_(self, RingElement right):656"""657Computes the quotient of ``self`` and ``right``.658659INPUT:660661- ``right`` - RingElement that is the divisor662663OUTPUT:664665- Quotient of ``self`` and ``right``666667EXAMPLES::668669sage: K.<x,y,z> = Frac(ZZ['x,y,z'])670sage: a = (x+1)*(x+y)/(z-3)671sage: b = (x+y)/(z-1)672sage: a/b673(x*z - x + z - 1)/(z - 3)674"""675snum = (<FractionFieldElement> right).__numerator676sden = (<FractionFieldElement> right).__denominator677678if snum.is_zero():679raise ZeroDivisionError, "fraction field element division by zero"680681rightinv = self.__class__(self._parent, sden, snum,682coerce=True, reduce=False)683684return self._mul_(rightinv)685686def __int__(self):687"""688EXAMPLES::689690sage: K = Frac(ZZ['x'])691sage: int(K(-3))692-3693sage: K.<x> = Frac(RR['x'])694sage: x/x695x/x696sage: int(x/x)6971698sage: int(K(.5))6990700"""701if self.__denominator != 1:702self.reduce()703if self.__denominator == 1:704return int(self.__numerator)705else:706raise TypeError, "denominator must equal 1"707708def _integer_(self, Z=ZZ):709"""710EXAMPLES::711712sage: K = Frac(ZZ['x'])713sage: K(5)._integer_()7145715sage: K.<x> = Frac(RR['x'])716sage: ZZ(2*x/x)7172718"""719if self.__denominator != 1:720self.reduce()721if self.__denominator == 1:722return Z(self.__numerator)723raise TypeError, "no way to coerce to an integer."724725def _rational_(self, Q=QQ):726"""727EXAMPLES::728729sage: K.<x> = Frac(QQ['x'])730sage: K(1/2)._rational_()7311/2732sage: K(1/2 + x/x)._rational_()7333/2734"""735return Q(self.__numerator) / Q(self.__denominator)736737def __long__(self):738"""739EXAMPLES::740741sage: K.<x> = Frac(QQ['x'])742sage: long(K(3))7433L744sage: long(K(3/5))7450L746"""747return long(int(self))748749def __pow__(self, right, dummy):750r"""751Returns self raised to the `right^{th}` power.752753Note that we need to check whether or not right is negative so we754don't set __numerator or __denominator to an element of the755fraction field instead of the underlying ring.756757EXAMPLES::758759sage: R = QQ['x','y']760sage: FR = R.fraction_field()761sage: x,y = FR.gens()762sage: a = x^2; a763x^2764sage: type(a.numerator())765<type 'sage.rings.polynomial.multi_polynomial_libsingular.MPolynomial_libsingular'>766sage: type(a.denominator())767<type 'sage.rings.polynomial.multi_polynomial_libsingular.MPolynomial_libsingular'>768sage: a = x^(-2); a7691/x^2770sage: type(a.numerator())771<type 'sage.rings.polynomial.multi_polynomial_libsingular.MPolynomial_libsingular'>772sage: type(a.denominator())773<type 'sage.rings.polynomial.multi_polynomial_libsingular.MPolynomial_libsingular'>774sage: x^07751776sage: ((x+y)/(x-y))^2777(x^2 + 2*x*y + y^2)/(x^2 - 2*x*y + y^2)778sage: ((x+y)/(x-y))^-2779(x^2 - 2*x*y + y^2)/(x^2 + 2*x*y + y^2)780sage: ((x+y)/(x-y))^07811782"""783snum = (<FractionFieldElement> self).__numerator784sden = (<FractionFieldElement> self).__denominator785if right == 0:786R = self.parent().ring()787return self.__class__(self.parent(),788R.one_element(), R.one_element(),789coerce=False, reduce=False)790elif right > 0:791return self.__class__(self.parent(),792snum**right, sden**right,793coerce=False, reduce=False)794else:795right = -right796return self.__class__(self.parent(),797sden**right, snum**right,798coerce=False, reduce=False)799800def __neg__(self):801"""802EXAMPLES::803804sage: K.<t> = Frac(GF(5)['t'])805sage: f = (t^2+t)/(t+2); f806(t^2 + t)/(t + 2)807sage: -f808(4*t^2 + 4*t)/(t + 2)809"""810return self.__class__(self._parent,811-self.__numerator, self.__denominator,812coerce=False, reduce=False)813814def __abs__(self):815"""816EXAMPLES::817818sage: from sage.rings.fraction_field_element import FractionFieldElement819sage: abs(FractionFieldElement(QQ, -2, 3, coerce=False))8202/3821"""822return abs(self.__numerator)/abs(self.__denominator)823824def __invert__(self):825"""826EXAMPLES::827828sage: K.<t> = Frac(GF(7)['t'])829sage: f = (t^2+5)/(t-1)830sage: ~f831(t + 6)/(t^2 + 5)832"""833if self.is_zero():834raise ZeroDivisionError, "Cannot invert 0"835return self.__class__(self._parent,836self.__denominator, self.__numerator, coerce=False, reduce=False)837838def __float__(self):839"""840EXAMPLES::841842sage: K.<x,y> = Frac(ZZ['x,y'])843sage: float(x/x + y/y)8442.0845"""846return float(self.__numerator) / float(self.__denominator)847848def __richcmp__(left, right, int op):849"""850EXAMPLES::851852sage: K.<x,y> = Frac(ZZ['x,y'])853sage: x > y854True855sage: 1 > y856False857"""858return (<Element>left)._richcmp(right, op)859860cdef int _cmp_c_impl(self, Element other) except -2:861"""862EXAMPLES::863864sage: K.<t> = Frac(GF(7)['t'])865sage: t/t == 1866True867sage: t+1/t == (t^2+1)/t868True869sage: t == t/5870False871"""872return cmp(self.__numerator * \873(<FractionFieldElement>other).__denominator,874self.__denominator*(<FractionFieldElement>other).__numerator)875876def valuation(self, v=None):877"""878Return the valuation of self, assuming that the numerator and879denominator have valuation functions defined on them.880881EXAMPLES::882883sage: x = PolynomialRing(RationalField(),'x').gen()884sage: f = (x^3 + x)/(x^2 - 2*x^3)885sage: f886(x^2 + 1)/(-2*x^2 + x)887sage: f.valuation()888-1889sage: f.valuation(x^2+1)8901891"""892return self.__numerator.valuation(v) - self.__denominator.valuation(v)893894def __nonzero__(self):895"""896Returns True if this element is nonzero.897898EXAMPLES::899900sage: F = ZZ['x,y'].fraction_field()901sage: x,y = F.gens()902sage: t = F(0)/x903sage: t.__nonzero__()904False905906::907908sage: (1/x).__nonzero__()909True910"""911return not self.__numerator.is_zero()912913def is_zero(self):914"""915Returns True if this element is equal to zero.916917EXAMPLES::918919sage: F = ZZ['x,y'].fraction_field()920sage: x,y = F.gens()921sage: t = F(0)/x922sage: t.is_zero()923True924sage: u = 1/x - 1/x925sage: u.is_zero()926True927sage: u.parent() is F928True929"""930return self.__numerator.is_zero()931932def is_one(self):933"""934Returns True if this element is equal to one.935936EXAMPLES::937938sage: F = ZZ['x,y'].fraction_field()939sage: x,y = F.gens()940sage: (x/x).is_one()941True942sage: (x/y).is_one()943False944"""945return self.__numerator == self.__denominator946947def _symbolic_(self, ring):948return ring(self.__numerator)/ring(self.__denominator)949950def __reduce__(self):951"""952EXAMPLES::953954sage: F = ZZ['x,y'].fraction_field()955sage: f = F.random_element()956sage: loads(f.dumps()) == f957True958"""959return (make_element,960(self._parent, self.__numerator, self.__denominator))961962963class FractionFieldElement_1poly_field(FractionFieldElement):964"""965A fraction field element where the parent is the fraction field of a univariate polynomial ring.966967Many of the functions here are included for coherence with number fields.968"""969def is_integral(self):970"""971Returns whether this element is actually a polynomial.972973EXAMPLES::974975sage: R.<t> = GF(5)[]976sage: K = R.fraction_field977"""978if self.__denominator != 1:979self.reduce()980return self.__denominator == 1981982def support(self):983"""984Returns a sorted list of primes dividing either the numerator or denominator of this element.985986EXAMPLES::987988sage: R.<t> = QQ[]989sage: 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)990sage: h.support()991[t - 1, t + 3, t^2 + 2, t^2 + t + 1, t^4 - 2]992"""993L = [fac[0] for fac in self.numerator().factor()] + [fac[0] for fac in self.denominator().factor()]994L.sort()995return L996997998999def make_element(parent, numerator, denominator):1000"""1001Used for unpickling FractionFieldElement objects (and subclasses).10021003EXAMPLES::10041005sage: from sage.rings.fraction_field_element import make_element1006sage: R = ZZ['x,y']1007sage: x,y = R.gens()1008sage: F = R.fraction_field()1009sage: make_element(F, 1+x, 1+y)1010(x + 1)/(y + 1)1011"""10121013return parent._element_class(parent, numerator, denominator)10141015def make_element_old(parent, cdict):1016"""1017Used for unpickling old FractionFieldElement pickles.10181019EXAMPLES::10201021sage: from sage.rings.fraction_field_element import make_element_old1022sage: R.<x,y> = ZZ[]1023sage: F = R.fraction_field()1024sage: make_element_old(F, {'_FractionFieldElement__numerator':x+y,'_FractionFieldElement__denominator':x-y})1025(x + y)/(x - y)1026"""1027return FractionFieldElement(parent,1028cdict['_FractionFieldElement__numerator'],1029cdict['_FractionFieldElement__denominator'],1030coerce=False, reduce=False)1031103210331034