Path: blob/master/src/sage/rings/fraction_field_element.pyx
8817 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 a :class`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 :trac:`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 :trac:`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"""103Initialize ``self``.104105EXAMPLES::106107sage: from sage.rings.fraction_field_element import FractionFieldElement108sage: K.<x> = Frac(ZZ['x'])109sage: FractionFieldElement(K, x, 4)110x/4111sage: FractionFieldElement(K, x, x, reduce=False)112x/x113sage: f = FractionFieldElement(K, 'hi', 1, coerce=False, reduce=False)114sage: f.numerator()115'hi'116117sage: x = var('x')118sage: K((x + 1)/(x^2 + x + 1))119(x + 1)/(x^2 + x + 1)120sage: K(355/113)121355/113122123"""124FieldElement.__init__(self, parent)125if coerce:126self.__numerator = parent.ring()(numerator)127self.__denominator = parent.ring()(denominator)128else:129self.__numerator = numerator130self.__denominator = denominator131if reduce and parent.is_exact():132try:133self.reduce()134except ArithmeticError:135pass136if self.__denominator.is_zero():137raise ZeroDivisionError, "fraction field element division by zero"138139def _im_gens_(self, codomain, im_gens):140"""141EXAMPLES::142143sage: F = ZZ['x,y'].fraction_field()144sage: x,y = F.gens()145sage: K = GF(7)['a,b'].fraction_field()146sage: a,b = K.gens()147148::149150sage: phi = F.hom([a+b, a*b], K)151sage: phi(x+y) # indirect doctest152a*b + a + b153154::155156sage: (x^2/y)._im_gens_(K, [a+b, a*b])157(a^2 + 2*a*b + b^2)/(a*b)158sage: (x^2/y)._im_gens_(K, [a, a*b])159a/b160"""161nnum = codomain.coerce(self.__numerator._im_gens_(codomain, im_gens))162nden = codomain.coerce(self.__denominator._im_gens_(codomain, im_gens))163return codomain.coerce(nnum/nden)164165def reduce(self):166"""167Divides out the gcd of the numerator and denominator.168169Automatically called for exact rings, but because it may be170numerically unstable for inexact rings it must be called manually171in that case.172173EXAMPLES::174175sage: R.<x> = RealField(10)[]176sage: f = (x^2+2*x+1)/(x+1); f177(x^2 + 2.0*x + 1.0)/(x + 1.0)178sage: f.reduce(); f179x + 1.0180"""181try:182g = self.__numerator.gcd(self.__denominator)183if not g.is_unit():184num, _ = self.__numerator.quo_rem(g)185den, _ = self.__denominator.quo_rem(g)186else:187num = self.__numerator188den = self.__denominator189if not den.is_one() and den.is_unit():190try:191num *= den.inverse_of_unit()192den = den.parent().one_element()193except Exception:194pass195self.__numerator = num196self.__denominator = den197except AttributeError:198raise ArithmeticError, "unable to reduce because lack of gcd or quo_rem algorithm"199except TypeError:200raise ArithmeticError, "unable to reduce because gcd algorithm doesn't work on input"201except NotImplementedError:202raise ArithmeticError("unable to reduce because gcd algorithm not implemented on input")203204def __copy__(self):205"""206Make a copy of ``self``.207208EXAMPLES::209210sage: R.<x,y> = ZZ[]211sage: f = x/y+1; f212(x + y)/y213sage: copy(f)214(x + y)/y215"""216return self.__class__(self._parent, self.__numerator,217self.__denominator, coerce=False, reduce=False)218219def numerator(self):220"""221Return the numerator of ``self``.222223EXAMPLES::224225sage: R.<x,y> = ZZ[]226sage: f = x/y+1; f227(x + y)/y228sage: f.numerator()229x + y230"""231return self.__numerator232233def denominator(self):234"""235Return the denominator of ``self``.236237EXAMPLES::238239sage: R.<x,y> = ZZ[]240sage: f = x/y+1; f241(x + y)/y242sage: f.denominator()243y244"""245return self.__denominator246247248def is_square(self,root=False):249"""250Returns whether or not ``self`` is a perfect square. If the optional251argument ``root`` is ``True``, then also returns a square root (or252``None``, if the fraction field element is not square).253254INPUT:255256- ``root`` -- whether or not to also return a square257root (default: ``False``)258259OUTPUT:260261- ``bool`` - whether or not a square262263- ``object`` - (optional) an actual square root if264found, and None otherwise.265266EXAMPLES::267268sage: R.<t> = QQ[]269sage: (1/t).is_square()270False271sage: (1/t^6).is_square()272True273sage: ((1+t)^4/t^6).is_square()274True275sage: (4*(1+t)^4/t^6).is_square()276True277sage: (2*(1+t)^4/t^6).is_square()278False279sage: ((1+t)/t^6).is_square()280False281282sage: (4*(1+t)^4/t^6).is_square(root=True)283(True, (2*t^2 + 4*t + 2)/t^3)284sage: (2*(1+t)^4/t^6).is_square(root=True)285(False, None)286287sage: R.<x> = QQ[]288sage: a = 2*(x+1)^2 / (2*(x-1)^2); a289(2*x^2 + 4*x + 2)/(2*x^2 - 4*x + 2)290sage: a.numerator().is_square()291False292sage: a.is_square()293True294sage: (0/x).is_square()295True296"""297a = self.numerator()298b = self.denominator()299if not root:300return (a*b).is_square( root = False )301is_sqr, sq_rt = (a*b).is_square( root = True )302if is_sqr:303return True, self._parent( sq_rt/b )304return False, None305306307def __hash__(self):308"""309This function hashes in a special way to ensure that generators of310a ring `R` and generators of a fraction field of `R` have the same311hash. This enables them to be used as keys interchangeably in a312dictionary (since ``==`` will claim them equal). This is particularly313useful for methods like ``subs`` on ``ParentWithGens`` if you are314passing a dictionary of substitutions.315316EXAMPLES::317318sage: R.<x>=ZZ[]319sage: hash(R.0)==hash(FractionField(R).0)320True321sage: ((x+1)/(x^2+1)).subs({x:1})3221323sage: d={x:1}324sage: d[FractionField(R).0]3251326sage: R.<x>=QQ[] # this probably has a separate implementation from ZZ[]327sage: hash(R.0)==hash(FractionField(R).0)328True329sage: d={x:1}330sage: d[FractionField(R).0]3311332sage: R.<x,y,z>=ZZ[] # this probably has a separate implementation from ZZ[]333sage: hash(R.0)==hash(FractionField(R).0)334True335sage: d={x:1}336sage: d[FractionField(R).0]3371338sage: R.<x,y,z>=QQ[] # this probably has a separate implementation from ZZ[]339sage: hash(R.0)==hash(FractionField(R).0)340True341sage: ((x+1)/(x^2+1)).subs({x:1})3421343sage: d={x:1}344sage: d[FractionField(R).0]3451346sage: hash(R(1)/R(2))==hash(1/2)347True348"""349# This is same algorithm as used for members of QQ350#cdef long n, d351n = hash(self.__numerator)352d = hash(self.__denominator)353if d == 1:354return n355n = n ^ d356if n == -1:357return -2358return n359360def __call__(self, *x, **kwds):361"""362Evaluate the fraction at the given arguments. This assumes that a363call function is defined for the numerator and denominator.364365EXAMPLES::366367sage: x = PolynomialRing(RationalField(),'x',3).gens()368sage: f = x[0] + x[1] - 2*x[1]*x[2]369sage: f370-2*x1*x2 + x0 + x1371sage: f(1,2,5)372-17373sage: h = f /(x[1] + x[2])374sage: h375(-2*x1*x2 + x0 + x1)/(x1 + x2)376sage: h(1,2,5)377-17/7378sage: h(x0=1)379(-2*x1*x2 + x1 + 1)/(x1 + x2)380"""381return self.__numerator(*x, **kwds) / self.__denominator(*x, **kwds)382383def _is_atomic(self):384"""385EXAMPLES::386387sage: K.<x> = Frac(ZZ['x'])388sage: x._is_atomic()389True390sage: f = 1/(x+1)391sage: f._is_atomic()392False393"""394return self.__numerator._is_atomic() and self.__denominator._is_atomic()395396def _repr_(self):397"""398Return a string representation of ``self``.399400EXAMPLES::401402sage: K.<x> = Frac(ZZ['x'])403sage: repr(x+1) # indirect doctest404'x + 1'405sage: repr((x+1)/(x-1))406'(x + 1)/(x - 1)'407sage: repr(1/(x-1))408'1/(x - 1)'409sage: repr(1/x)410'1/x'411"""412if self.is_zero():413return "0"414s = "%s"%self.__numerator415if self.__denominator != 1:416denom_string = str( self.__denominator )417if self.__denominator._is_atomic() and not ('*' in denom_string or '/' in denom_string):418s = "%s/%s"%(self.__numerator._coeff_repr(no_space=False),denom_string)419else:420s = "%s/(%s)"%(self.__numerator._coeff_repr(no_space=False),denom_string)421return s422423def _latex_(self):424r"""425Return a latex representation of this fraction field element.426427EXAMPLES::428429sage: R = PolynomialRing(QQ, 'x')430sage: F = R.fraction_field()431sage: x = F.gen()432sage: a = x^2 / 1433sage: latex(a) # indirect doctest434x^{2}435sage: latex(x^2/(x^2+1))436\frac{x^{2}}{x^{2} + 1}437sage: a = 1/x438sage: latex(a)439\frac{1}{x}440441TESTS::442443sage: R = RR['x'] # Inexact, so no reduction.444sage: F = Frac(R)445sage: from sage.rings.fraction_field_element import FractionFieldElement446sage: z = FractionFieldElement(F, 0, R.gen(), coerce=False)447sage: z.numerator() == 0448True449sage: z.denominator() == R.gen()450True451sage: latex(z) # indirect doctest4520453"""454if self.is_zero():455return "0"456if self.__denominator == 1:457return latex.latex(self.__numerator)458return "\\frac{%s}{%s}"%(latex.latex(self.__numerator),459latex.latex(self.__denominator))460461def _magma_init_(self, magma):462"""463Return a string representation of ``self`` Magma can understand.464465EXAMPLES::466467sage: R.<x> = ZZ[]468sage: magma((x^2 + x + 1)/(x + 1)) # optional - magma # indirect doctest469(x^2 + x + 1)/(x + 1)470471::472473sage: R.<x,y> = QQ[]474sage: magma((x+y)/x) # optional - magma475(x + y)/x476"""477pgens = magma(self._parent).gens()478479s = self._repr_()480for i, j in zip(self._parent.variable_names(), pgens):481s = s.replace(i, j.name())482483return s484485cpdef ModuleElement _add_(self, ModuleElement right):486"""487Computes the sum of ``self`` and ``right``.488489INPUT:490491- ``right`` -- ``ModuleElement`` to add to ``self``492493OUTPUT:494495- Sum of ``self`` and ``right``496497EXAMPLES::498499sage: K.<x,y> = Frac(ZZ['x,y'])500sage: x+y # indirect doctest501x + y502sage: 1/x + 1/y503(x + y)/(x*y)504sage: 1/x + 1/(x*y)505(y + 1)/(x*y)506sage: Frac(CDF['x']).gen() + 3507x + 3.0508"""509rnum = self.__numerator510rden = self.__denominator511snum = (<FractionFieldElement> right).__numerator512sden = (<FractionFieldElement> right).__denominator513514if (rnum.is_zero()):515return <FractionFieldElement> right516if (snum.is_zero()):517return self518519if self._parent.is_exact():520try:521d = rden.gcd(sden)522if d.is_unit():523return self.__class__(self._parent, rnum*sden + rden*snum,524rden*sden, coerce=False, reduce=False)525else:526rden = rden // d527sden = sden // d528tnum = rnum * sden + rden * snum529if tnum.is_zero():530return self.__class__(self._parent, tnum,531self._parent.ring().one_element(), coerce=False,532reduce=False)533else:534tden = self.__denominator * sden535e = tnum.gcd(d)536if not e.is_unit():537tnum = tnum // e538tden = tden // e539if not tden.is_one() and tden.is_unit():540try:541tnum = tnum * tden.inverse_of_unit()542tden = self._parent.ring().one_element()543except AttributeError:544pass545except NotImplementedError:546pass547return self.__class__(self._parent, tnum, tden,548coerce=False, reduce=False)549except AttributeError:550pass551except NotImplementedError:552pass553except TypeError:554pass555556rnum = self.__numerator557rden = self.__denominator558snum = (<FractionFieldElement> right).__numerator559sden = (<FractionFieldElement> right).__denominator560561return self.__class__(self._parent, rnum*sden + rden*snum, rden*sden,562coerce=False, reduce=False)563564cpdef ModuleElement _sub_(self, ModuleElement right):565"""566Computes the difference of ``self`` and ``right``.567568INPUT:569570- ``right`` - ``ModuleElement`` to subtract from ``self``571572OUTPUT:573574- Difference of ``self`` and ``right``575576EXAMPLES::577578sage: K.<t> = Frac(GF(7)['t'])579sage: t - 1/t # indirect doctest580(t^2 + 6)/t581"""582return self._add_(-right)583584cpdef RingElement _mul_(self, RingElement right):585"""586Computes the product of ``self`` and ``right``.587588INPUT:589590- ``right`` - ``RingElement`` to multiply with ``self``591592OUTPUT:593594- Product of ``self`` and ``right``595596EXAMPLES::597598sage: K.<t> = Frac(GF(7)['t'])599sage: a = t/(1+t)600sage: b = 3/t601sage: a*b # indirect doctest6023/(t + 1)603"""604rnum = self.__numerator605rden = self.__denominator606snum = (<FractionFieldElement> right).__numerator607sden = (<FractionFieldElement> right).__denominator608609if (rnum.is_zero() or snum.is_zero()):610return self._parent.zero_element()611612if self._parent.is_exact():613try:614d1 = rnum.gcd(sden)615d2 = snum.gcd(rden)616if not d1.is_unit():617rnum = rnum // d1618sden = sden // d1619if not d2.is_unit():620rden = rden // d2621snum = snum // d2622tnum = rnum * snum623tden = rden * sden624if not tden.is_one() and tden.is_unit():625try:626tnum = tnum * tden.inverse_of_unit()627tden = self._parent.ring().one_element()628except AttributeError:629pass630except NotImplementedError:631pass632return self.__class__(self._parent, tnum, tden,633coerce=False, reduce=False)634except AttributeError:635pass636except NotImplementedError:637pass638except TypeError:639pass640641rnum = self.__numerator642rden = self.__denominator643snum = (<FractionFieldElement> right).__numerator644sden = (<FractionFieldElement> right).__denominator645646return self.__class__(self._parent, rnum * snum, rden * sden,647coerce=False, reduce=False)648649cpdef RingElement _div_(self, RingElement right):650"""651Computes the quotient of ``self`` and ``right``.652653INPUT:654655- ``right`` -- ``RingElement`` that is the divisor656657OUTPUT:658659Quotient of ``self`` and ``right``660661EXAMPLES::662663sage: K.<x,y,z> = Frac(ZZ['x,y,z'])664sage: a = (x+1)*(x+y)/(z-3) # indirect doctest665sage: b = (x+y)/(z-1)666sage: a/b667(x*z - x + z - 1)/(z - 3)668"""669snum = (<FractionFieldElement> right).__numerator670sden = (<FractionFieldElement> right).__denominator671672if snum.is_zero():673raise ZeroDivisionError, "fraction field element division by zero"674675rightinv = self.__class__(self._parent, sden, snum,676coerce=True, reduce=False)677678return self._mul_(rightinv)679680def __int__(self):681"""682EXAMPLES::683684sage: K = Frac(ZZ['x'])685sage: int(K(-3))686-3687sage: K.<x> = Frac(RR['x'])688sage: x/x689x/x690sage: int(x/x)6911692sage: int(K(.5))6930694"""695if self.__denominator != 1:696self.reduce()697if self.__denominator == 1:698return int(self.__numerator)699else:700raise TypeError, "denominator must equal 1"701702def _integer_(self, Z=ZZ):703"""704EXAMPLES::705706sage: K = Frac(ZZ['x'])707sage: K(5)._integer_()7085709sage: K.<x> = Frac(RR['x'])710sage: ZZ(2*x/x)7112712"""713if self.__denominator != 1:714self.reduce()715if self.__denominator == 1:716return Z(self.__numerator)717raise TypeError, "no way to coerce to an integer."718719def _rational_(self, Q=QQ):720"""721EXAMPLES::722723sage: K.<x> = Frac(QQ['x'])724sage: K(1/2)._rational_()7251/2726sage: K(1/2 + x/x)._rational_()7273/2728"""729return Q(self.__numerator) / Q(self.__denominator)730731def __long__(self):732"""733EXAMPLES::734735sage: K.<x> = Frac(QQ['x'])736sage: long(K(3))7373L738sage: long(K(3/5))7390L740"""741return long(int(self))742743def __pow__(self, right, dummy):744r"""745Returns self raised to the `right^{th}` power.746747Note that we need to check whether or not right is negative so we748don't set ``__numerator`` or ``__denominator`` to an element of the749fraction field instead of the underlying ring.750751EXAMPLES::752753sage: R = QQ['x','y']754sage: FR = R.fraction_field()755sage: x,y = FR.gens()756sage: a = x^2; a757x^2758sage: type(a.numerator())759<type 'sage.rings.polynomial.multi_polynomial_libsingular.MPolynomial_libsingular'>760sage: type(a.denominator())761<type 'sage.rings.polynomial.multi_polynomial_libsingular.MPolynomial_libsingular'>762sage: a = x^(-2); a7631/x^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: x^07691770sage: ((x+y)/(x-y))^2771(x^2 + 2*x*y + y^2)/(x^2 - 2*x*y + y^2)772sage: ((x+y)/(x-y))^-2773(x^2 - 2*x*y + y^2)/(x^2 + 2*x*y + y^2)774sage: ((x+y)/(x-y))^07751776"""777snum = (<FractionFieldElement> self).__numerator778sden = (<FractionFieldElement> self).__denominator779if right == 0:780R = self.parent().ring()781return self.__class__(self.parent(),782R.one_element(), R.one_element(),783coerce=False, reduce=False)784elif right > 0:785return self.__class__(self.parent(),786snum**right, sden**right,787coerce=False, reduce=False)788else:789right = -right790return self.__class__(self.parent(),791sden**right, snum**right,792coerce=False, reduce=False)793794def __neg__(self):795"""796EXAMPLES::797798sage: K.<t> = Frac(GF(5)['t'])799sage: f = (t^2+t)/(t+2); f800(t^2 + t)/(t + 2)801sage: -f802(4*t^2 + 4*t)/(t + 2)803"""804return self.__class__(self._parent,805-self.__numerator, self.__denominator,806coerce=False, reduce=False)807808def __abs__(self):809"""810EXAMPLES::811812sage: from sage.rings.fraction_field_element import FractionFieldElement813sage: abs(FractionFieldElement(QQ, -2, 3, coerce=False))8142/3815"""816return abs(self.__numerator)/abs(self.__denominator)817818def __invert__(self):819"""820EXAMPLES::821822sage: K.<t> = Frac(GF(7)['t'])823sage: f = (t^2+5)/(t-1)824sage: ~f825(t + 6)/(t^2 + 5)826"""827if self.is_zero():828raise ZeroDivisionError, "Cannot invert 0"829return self.__class__(self._parent,830self.__denominator, self.__numerator, coerce=False, reduce=False)831832def __float__(self):833"""834EXAMPLES::835836sage: K.<x,y> = Frac(ZZ['x,y'])837sage: float(x/x + y/y)8382.0839"""840return float(self.__numerator) / float(self.__denominator)841842def __richcmp__(left, right, int op):843"""844EXAMPLES::845846sage: K.<x,y> = Frac(ZZ['x,y'])847sage: x > y848True849sage: 1 > y850False851"""852return (<Element>left)._richcmp(right, op)853854cdef int _cmp_c_impl(self, Element other) except -2:855"""856EXAMPLES::857858sage: K.<t> = Frac(GF(7)['t'])859sage: t/t == 1860True861sage: t+1/t == (t^2+1)/t862True863sage: t == t/5864False865"""866return cmp(self.__numerator * \867(<FractionFieldElement>other).__denominator,868self.__denominator*(<FractionFieldElement>other).__numerator)869870def valuation(self, v=None):871"""872Return the valuation of ``self``, assuming that the numerator and873denominator have valuation functions defined on them.874875EXAMPLES::876877sage: x = PolynomialRing(RationalField(),'x').gen()878sage: f = (x^3 + x)/(x^2 - 2*x^3)879sage: f880(x^2 + 1)/(-2*x^2 + x)881sage: f.valuation()882-1883sage: f.valuation(x^2+1)8841885"""886return self.__numerator.valuation(v) - self.__denominator.valuation(v)887888def __nonzero__(self):889"""890Return ``True`` if this element is nonzero.891892EXAMPLES::893894sage: F = ZZ['x,y'].fraction_field()895sage: x,y = F.gens()896sage: t = F(0)/x897sage: t.__nonzero__()898False899900::901902sage: (1/x).__nonzero__()903True904"""905return not self.__numerator.is_zero()906907def is_zero(self):908"""909Return ``True`` if this element is equal to zero.910911EXAMPLES::912913sage: F = ZZ['x,y'].fraction_field()914sage: x,y = F.gens()915sage: t = F(0)/x916sage: t.is_zero()917True918sage: u = 1/x - 1/x919sage: u.is_zero()920True921sage: u.parent() is F922True923"""924return self.__numerator.is_zero()925926def is_one(self):927"""928Return ``True`` if this element is equal to one.929930EXAMPLES::931932sage: F = ZZ['x,y'].fraction_field()933sage: x,y = F.gens()934sage: (x/x).is_one()935True936sage: (x/y).is_one()937False938"""939return self.__numerator == self.__denominator940941def _symbolic_(self, ring):942"""943Return ``self`` as a fraction in the ring ``ring``. Used for944:func:`symbolic_expression` in creating a symbolic expression of945``self``.946947EXAMPLES::948949sage: F = ZZ['x,y'].fraction_field()950sage: x,y = F.gens()951sage: elt = (2*x + 2*y) / (3*x - 3*y); elt952(2*x + 2*y)/(3*x - 3*y)953sage: elt._symbolic_(SR)9542/3*(x + y)/(x - y)955sage: symbolic_expression(elt)9562/3*(x + y)/(x - y)957"""958return ring(self.__numerator)/ring(self.__denominator)959960def __reduce__(self):961"""962For pickling.963964EXAMPLES::965966sage: F = ZZ['x,y'].fraction_field()967sage: f = F.random_element()968sage: loads(f.dumps()) == f969True970"""971return (make_element,972(self._parent, self.__numerator, self.__denominator))973974975class FractionFieldElement_1poly_field(FractionFieldElement):976"""977A fraction field element where the parent is the fraction field of a978univariate polynomial ring.979980Many of the functions here are included for coherence with number fields.981"""982def is_integral(self):983"""984Returns whether this element is actually a polynomial.985986EXAMPLES::987988sage: R.<t> = QQ[]989sage: elt = (t^2 + t - 2) / (t + 2); elt # == (t + 2)*(t - 1)/(t + 2)990t - 1991sage: elt.is_integral()992True993sage: elt = (t^2 - t) / (t+2); elt # == t*(t - 1)/(t + 2)994(t^2 - t)/(t + 2)995sage: elt.is_integral()996False997"""998if self.denominator() != 1:999self.reduce()1000return self.denominator() == 110011002def support(self):1003"""1004Returns a sorted list of primes dividing either the numerator or1005denominator of this element.10061007EXAMPLES::10081009sage: R.<t> = QQ[]1010sage: 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)1011sage: h.support()1012[t - 1, t + 3, t^2 + 2, t^2 + t + 1, t^4 - 2]1013"""1014L = [fac[0] for fac in self.numerator().factor()] + [fac[0] for fac in self.denominator().factor()]1015L.sort()1016return L101710181019def make_element(parent, numerator, denominator):1020"""1021Used for unpickling :class:`FractionFieldElement` objects (and subclasses).10221023EXAMPLES::10241025sage: from sage.rings.fraction_field_element import make_element1026sage: R = ZZ['x,y']1027sage: x,y = R.gens()1028sage: F = R.fraction_field()1029sage: make_element(F, 1+x, 1+y)1030(x + 1)/(y + 1)1031"""10321033return parent._element_class(parent, numerator, denominator)10341035def make_element_old(parent, cdict):1036"""1037Used for unpickling old :class:`FractionFieldElement` pickles.10381039EXAMPLES::10401041sage: from sage.rings.fraction_field_element import make_element_old1042sage: R.<x,y> = ZZ[]1043sage: F = R.fraction_field()1044sage: make_element_old(F, {'_FractionFieldElement__numerator':x+y,'_FractionFieldElement__denominator':x-y})1045(x + y)/(x - y)1046"""1047return FractionFieldElement(parent,1048cdict['_FractionFieldElement__numerator'],1049cdict['_FractionFieldElement__denominator'],1050coerce=False, reduce=False)1051105210531054