Skip to content

Commit 150e48f

Browse files
Release Managervbraun
Release Manager
authored andcommitted
Trac #13629: provide xgcd for new polynomial rings through _xgcd_univariate_polynomial
Currently, to add xgcd functionality for a new polynomial ring, one needs to add a specialized subclass of !PolynomialElement. The attached patch allows rings to provide a {{{_xgcd_univariate_polynomial}}} method which will be called by !PolynomialElement to compute xgcds. This is similar to #10635, #13442. URL: http://trac.sagemath.org/13629 Reported by: saraedum Ticket author(s): Julian Rueth Reviewer(s): Peter Bruin, Bruno Grenet
2 parents 005bf73 + 828b5fc commit 150e48f

File tree

5 files changed

+211
-81
lines changed

5 files changed

+211
-81
lines changed

src/doc/en/thematic_tutorials/coercion_and_categories.rst

+1
Original file line numberDiff line numberDiff line change
@@ -132,6 +132,7 @@ This base class provides a lot more methods than a general parent::
132132
'_pseudo_fraction_field',
133133
'_random_nonzero_element',
134134
'_unit_ideal',
135+
'_xgcd_univariate_polynomial',
135136
'_zero_element',
136137
'_zero_ideal',
137138
'algebraic_closure',

src/sage/categories/fields.py

+83-1
Original file line numberDiff line numberDiff line change
@@ -207,7 +207,7 @@ def _gcd_univariate_polynomial(self, f, g):
207207
208208
INPUT:
209209
210-
- ``f``, ``g`` -- two polynomials defined over ``self``
210+
- ``f``, ``g`` -- two polynomials defined over ``self``
211211
212212
.. NOTE::
213213
@@ -368,6 +368,88 @@ def __pow__(self, n):
368368
from sage.modules.all import FreeModule
369369
return FreeModule(self, n)
370370

371+
def _xgcd_univariate_polynomial(self, left, right):
372+
r"""
373+
Return an extended gcd of ``left`` and ``right``.
374+
375+
INPUT:
376+
377+
- ``left``, ``right`` -- two polynomials over this field
378+
379+
OUTPUT:
380+
381+
Polynomials ``g``, ``u``, and ``v`` such that ``g`` is a
382+
greatest common divisor of ``left and ``right``, and such
383+
that ``g = u*left + v*right`` holds.
384+
385+
.. NOTE::
386+
387+
This is a helper method for
388+
:meth:`sage.rings.polynomial.polynomial_element.Polynomial.xgcd`.
389+
390+
EXAMPLES::
391+
392+
sage: P.<x> = QQ[]
393+
sage: F = (x^2 + 2)*x^3; G = (x^2+2)*(x-3)
394+
sage: g, u, v = QQ._xgcd_univariate_polynomial(F,G)
395+
sage: g, u, v
396+
(x^2 + 2, 1/27, -1/27*x^2 - 1/9*x - 1/3)
397+
sage: u*F + v*G
398+
x^2 + 2
399+
400+
::
401+
402+
sage: g, u, v = QQ._xgcd_univariate_polynomial(x,P(0)); g, u, v
403+
(x, 1, 0)
404+
sage: g == u*x + v*P(0)
405+
True
406+
sage: g, u, v = QQ._xgcd_univariate_polynomial(P(0),x); g, u, v
407+
(x, 0, 1)
408+
sage: g == u*P(0) + v*x
409+
True
410+
411+
TESTS:
412+
413+
We check that the behavior of xgcd with zero elements is
414+
compatible with gcd (:trac:`17671`)::
415+
416+
sage: R.<x> = QQbar[]
417+
sage: zero = R.zero()
418+
sage: zero.xgcd(2*x)
419+
(x, 0, 1/2)
420+
sage: (2*x).xgcd(zero)
421+
(x, 1/2, 0)
422+
sage: zero.xgcd(zero)
423+
(0, 0, 0)
424+
"""
425+
R = left.parent()
426+
zero = R.zero()
427+
one = R.one()
428+
if right.is_zero():
429+
if left.is_zero():
430+
return (zero, zero, zero)
431+
else:
432+
c = left.leading_coefficient()
433+
return (left/c, one/c, zero)
434+
elif left.is_zero():
435+
c = right.leading_coefficient()
436+
return (right/c, zero, one/c)
437+
438+
# Algorithm 3.2.2 of Cohen, GTM 138
439+
A = left
440+
B = right
441+
U = one
442+
G = A
443+
V1 = zero
444+
V3 = B
445+
while not V3.is_zero():
446+
Q, R = G.quo_rem(V3)
447+
G, U, V1, V3 = V3, V1, U-V1*Q, R
448+
V = (G-A*U)//B
449+
lc = G.leading_coefficient()
450+
return G/lc, U/lc, V/lc
451+
452+
371453
class ElementMethods:
372454
def euclidean_degree(self):
373455
r"""

src/sage/rings/polynomial/polynomial_element.pyx

+65-5
Original file line numberDiff line numberDiff line change
@@ -3774,17 +3774,17 @@ cdef class Polynomial(CommutativeAlgebraElement):
37743774
@coerce_binop
37753775
def gcd(self, other):
37763776
"""
3777-
Compute a greatest common divisor of this polynomial and ``other``.
3777+
Return a greatest common divisor of this polynomial and ``other``.
37783778
37793779
INPUT:
37803780
3781-
- ``other`` -- a polynomial in the same ring as this polynomial
3781+
- ``other`` -- a polynomial in the same ring as this polynomial
37823782
37833783
OUTPUT:
37843784
3785-
A greatest common divisor as a polynomial in the same ring as this
3786-
polynomial. Over a field, the return value will be a monic
3787-
polynomial.
3785+
A greatest common divisor as a polynomial in the same ring as
3786+
this polynomial. If the base ring is a field, the return value
3787+
is a monic polynomial.
37883788
37893789
.. NOTE::
37903790
@@ -6207,6 +6207,66 @@ cdef class Polynomial(CommutativeAlgebraElement):
62076207
"""
62086208
return self.parent().variable_name()
62096209

6210+
@coerce_binop
6211+
def xgcd(self, other):
6212+
r"""
6213+
Return an extended gcd of this polynomial and ``other``.
6214+
6215+
INPUT:
6216+
6217+
- ``other`` -- a polynomial in the same ring as this polynomial
6218+
6219+
OUTPUT:
6220+
6221+
A tuple ``(r, s, t)`` where ``r`` is a greatest common divisor
6222+
of this polynomial and ``other``, and ``s`` and ``t`` are such
6223+
that ``r = s*self + t*other`` holds.
6224+
6225+
.. NOTE::
6226+
6227+
The actual algorithm for computing the extended gcd depends on the
6228+
base ring underlying the polynomial ring. If the base ring defines
6229+
a method ``_xgcd_univariate_polynomial``, then this method will be
6230+
called (see examples below).
6231+
6232+
EXAMPLES::
6233+
6234+
sage: R.<x> = QQbar[]
6235+
sage: (2*x^2).gcd(2*x)
6236+
x
6237+
sage: R.zero().gcd(0)
6238+
0
6239+
sage: (2*x).gcd(0)
6240+
x
6241+
6242+
One can easily add xgcd functionality to new rings by providing a
6243+
method ``_xgcd_univariate_polynomial``::
6244+
6245+
sage: R.<x> = QQ[]
6246+
sage: S.<y> = R[]
6247+
sage: h1 = y*x
6248+
sage: h2 = y^2*x^2
6249+
sage: h1.xgcd(h2)
6250+
Traceback (most recent call last):
6251+
...
6252+
NotImplementedError: Univariate Polynomial Ring in x over Rational Field does not provide an xgcd implementation for univariate polynomials
6253+
sage: T.<x,y> = QQ[]
6254+
sage: def poor_xgcd(f,g):
6255+
....: ret = S(T(f).gcd(g))
6256+
....: if ret == f: return ret,S.one(),S.zero()
6257+
....: if ret == g: return ret,S.zero(),S.one()
6258+
....: raise NotImplementedError
6259+
sage: R._xgcd_univariate_polynomial = poor_xgcd
6260+
sage: h1.xgcd(h2)
6261+
(x*y, 1, 0)
6262+
sage: del R._xgcd_univariate_polynomial
6263+
6264+
"""
6265+
if hasattr(self.base_ring(), '_xgcd_univariate_polynomial'):
6266+
return self.base_ring()._xgcd_univariate_polynomial(self, other)
6267+
else:
6268+
raise NotImplementedError("%s does not provide an xgcd implementation for univariate polynomials"%self.base_ring())
6269+
62106270
def variables(self):
62116271
"""
62126272
Returns the tuple of variables occurring in this polynomial.

src/sage/rings/polynomial/polynomial_element_generic.py

-75
Original file line numberDiff line numberDiff line change
@@ -726,81 +726,6 @@ def quo_rem(self, other):
726726
R = R[:R.degree()] - (aaa*B[:B.degree()]).shift(diff_deg)
727727
return (Q, R)
728728

729-
@coerce_binop
730-
def xgcd(self, other):
731-
r"""
732-
Extended gcd of ``self`` and polynomial ``other``.
733-
734-
INPUT:
735-
736-
- ``other`` -- a polynomial defined over the same ring as ``self``
737-
738-
OUTPUT:
739-
740-
Polynomials ``g``, ``u``, and ``v`` such that ``g = u * self + v * other``.
741-
742-
EXAMPLES::
743-
744-
sage: P.<x> = QQ[]
745-
sage: F = (x^2 + 2)*x^3; G = (x^2+2)*(x-3)
746-
sage: g, u, v = F.xgcd(G)
747-
sage: g, u, v
748-
(x^2 + 2, 1/27, -1/27*x^2 - 1/9*x - 1/3)
749-
sage: u*F + v*G
750-
x^2 + 2
751-
752-
::
753-
754-
sage: g, u, v = x.xgcd(P(0)); g, u, v
755-
(x, 1, 0)
756-
sage: g == u*x + v*P(0)
757-
True
758-
sage: g, u, v = P(0).xgcd(x); g, u, v
759-
(x, 0, 1)
760-
sage: g == u*P(0) + v*x
761-
True
762-
763-
TESTS:
764-
765-
We check that the behavior of xgcd with zero elements is compatible with
766-
gcd (:trac:`17671`)::
767-
768-
sage: R.<x> = QQbar[]
769-
sage: zero = R.zero()
770-
sage: zero.xgcd(2*x)
771-
(x, 0, 1/2)
772-
sage: (2*x).xgcd(zero)
773-
(x, 1/2, 0)
774-
sage: zero.xgcd(zero)
775-
(0, 0, 0)
776-
"""
777-
R = self.parent()
778-
zero = R.zero()
779-
one = R.one()
780-
if other.is_zero():
781-
if self.is_zero():
782-
return (zero, zero, zero)
783-
else:
784-
c = self.leading_coefficient()
785-
return (self/c, one/c, zero)
786-
elif self.is_zero():
787-
c = other.leading_coefficient()
788-
return (other/c, zero, one/c)
789-
790-
# Algorithm 3.2.2 of Cohen, GTM 138
791-
A = self
792-
B = other
793-
U = one
794-
G = A
795-
V1 = zero
796-
V3 = B
797-
while not V3.is_zero():
798-
Q, R = G.quo_rem(V3)
799-
G, U, V1, V3 = V3, V1, U-V1*Q, R
800-
V = (G-A*U)//B
801-
lc = G.leading_coefficient()
802-
return G/lc, U/lc, V/lc
803-
804729

805730
class Polynomial_generic_sparse_field(Polynomial_generic_sparse, Polynomial_generic_field):
806731
"""

src/sage/rings/ring.pyx

+62
Original file line numberDiff line numberDiff line change
@@ -2228,6 +2228,68 @@ cdef class Field(PrincipalIdealDomain):
22282228
a = a.monic()
22292229
return a
22302230

2231+
def _xgcd_univariate_polynomial(self, a, b):
2232+
"""
2233+
Return an extended gcd of ``a`` and ``b``.
2234+
2235+
INPUT:
2236+
2237+
- ``a``, ``b`` -- two univariate polynomials
2238+
2239+
OUTPUT:
2240+
2241+
A tuple ``(d, u, v)`` of polynomials such that ``d`` is the
2242+
greatest common divisor (monic or zero) of ``a`` and ``b``,
2243+
and ``u``, ``v`` satisfy ``d = u*a + v*b``.
2244+
2245+
.. WARNING:
2246+
2247+
If the base ring is inexact, the results may not be
2248+
entirely stable.
2249+
2250+
ALGORITHM:
2251+
2252+
This uses the extended Euclidean algorithm; see for example
2253+
[Cohen]_, Algorithm 3.2.2.
2254+
2255+
REFERENCES:
2256+
2257+
.. [Cohen] H. Cohen, A Course in Computational Algebraic
2258+
Number Theory. Graduate Texts in Mathematics 138.
2259+
Springer-Verlag, 1996.
2260+
2261+
TESTS::
2262+
2263+
sage: for A in (RR, CC, QQbar):
2264+
....: g = A._xgcd_univariate_polynomial
2265+
....: R.<x> = A[]
2266+
....: z, h = R(0), R(1/2)
2267+
....: assert(g(2*x, 2*x^2) == (x, h, z) and
2268+
....: g(z, 2*x) == (x, z, h) and
2269+
....: g(2*x, z) == (x, h, z) and
2270+
....: g(z, z) == (z, z, z))
2271+
2272+
"""
2273+
R = a.parent()
2274+
zero = R.zero()
2275+
if not b:
2276+
if not a:
2277+
return (zero, zero, zero)
2278+
c = ~a.leading_coefficient()
2279+
return (c*a, R(c), zero)
2280+
elif not a:
2281+
c = ~b.leading_coefficient()
2282+
return (c*b, zero, R(c))
2283+
(u, d, v1, v3) = (R.one(), a, zero, b)
2284+
while v3:
2285+
q, r = d.quo_rem(v3)
2286+
(u, d, v1, v3) = (v1, v3, u - v1*q, r)
2287+
v = (d - a*u) // b
2288+
if d:
2289+
c = ~d.leading_coefficient()
2290+
d, u, v = c*d, c*u, c*v
2291+
return d, u, v
2292+
22312293

22322294
cdef class Algebra(Ring):
22332295
"""

0 commit comments

Comments
 (0)