Skip to content

Commit c1f3652

Browse files
author
Release Manager
committed
gh-37443: Fix bugs and regression in performance of `any_root()` This pull request aims to fix some issues which came with the refactoring of `any_root()` in #37170. I am afraid that apart from a little more cleaning up and verbose comments, this "fix" really relies on `f.roots()` rather than `f.any_root()` when roots over extension fields are computed. The core problem currently is that the proper workflow should be the following: - Given a polynomial $f$, find an irreducible factor of smallest degree. When `ring` is `None` and `degree` is `None`, only linear factors are searched for. This code is unchanged and works fast. - When `degree` is `None` and `ring` is not `None`, the old code used to perform `f.change_ring(ring).any_root()` and look for a linear factor in `ring`. The issue is when `ring` is a field with a large degree, arithmetic in this ring is very slow and root finding is very slow. - When `degree` is not `None` then an irreducible of degree `degree` is searched for, $f$, and then a root is found in an extension ring. This is either the user supplied `ring`, or it is in the extension `self.base_ring().extension(degree)`. Again, as this extension could be large, this can be slow. Now, the reason that there's a regression in speed, as explained in #37359, is that the method before refactoring simply called `f.roots(ring)[0]` when working for these extensions. As the root finding is always performed by a C library this is much faster than the python code which we have for `any_root()` and so the more "pure" method I wrote simply is slower. The real issue holding this method back is that if we are in some field $GF(p^k)$ and `ring` is some extension $GF(p^{nk})$ then we are not guaranteed a coercion between these rings with the current coercion system. Furthermore, if we extend the base to some $GF(p^{dk})$ with $dk$ dividing $nk$ we also have no certainty of a coercion from this minimal extension to the user supplied ring. As a result, a "good" fix is delayed until we figure out finite field coercion. **Issue to come**. Because of all of this, this PR goes back to the old behaviour, while keeping the refactored and easier to read code. I hope one day in the future to fix this. Fixes #37359 Fixes #37442 Fixes #37445 Fixes #37471 **EDIT** I have labeled this as critical as the slow down in the most extreme cases causes code to complete hang when looking for a root. See for example: #37442. I do not want #37170 to be merged into 10.3 without this patch. Additionally, a typo meant that a bug was introduced in the CZ splitting. This has also been fixed. Also, also in the refactoring the function behaved as intended (rather than the strange behaviour from before) which exposed an issue #37471 which I have fixed. URL: #37443 Reported by: Giacomo Pope Reviewer(s): Giacomo Pope, grhkm21, Lorenz Panny, Travis Scrimshaw
2 parents 65a64ea + b8da002 commit c1f3652

File tree

4 files changed

+139
-51
lines changed

4 files changed

+139
-51
lines changed

src/sage/rings/finite_rings/conway_polynomials.py

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -239,7 +239,7 @@ def polynomial(self, n):
239239
# TODO: something like the following
240240
# gcds = [n.gcd(d) for d in self.nodes.keys()]
241241
# xi = { m: (...) for m in gcds }
242-
xi = {q: self.polynomial(n//q).any_root(K, n//q, assume_squarefree=True, assume_distinct_deg=True)
242+
xi = {q: self.polynomial(n//q).any_root(K, n//q, assume_squarefree=True, assume_equal_deg=True)
243243
for q in n.prime_divisors()}
244244

245245
# The following is needed to ensure that in the concrete instantiation

src/sage/rings/polynomial/polynomial_element.pyx

Lines changed: 126 additions & 46 deletions
Original file line numberDiff line numberDiff line change
@@ -134,7 +134,6 @@ from sage.categories.morphism cimport Morphism
134134
from sage.misc.superseded import deprecation_cython as deprecation, deprecated_function_alias
135135
from sage.misc.cachefunc import cached_method
136136

137-
138137
cpdef is_Polynomial(f) noexcept:
139138
"""
140139
Return ``True`` if ``f`` is of type univariate polynomial.
@@ -2151,6 +2150,18 @@ cdef class Polynomial(CommutativePolynomial):
21512150
True
21522151
sage: f % factor == 0
21532152
True
2153+
2154+
TESTS:
2155+
2156+
Ensure that :issue:`37445` is fixed::
2157+
2158+
sage: R.<x> = GF(13)[]
2159+
sage: def irr(d, R): return f.monic() if (f := R.random_element(d)).is_irreducible() else irr(d, R)
2160+
sage: f = prod(irr(6, R) for _ in range(10))
2161+
sage: irr = f._cantor_zassenhaus_split_to_irreducible(6)
2162+
sage: assert irr.degree() == 6
2163+
sage: assert f % irr == 0
2164+
sage: assert irr.is_irreducible()
21542165
"""
21552166
R = self.parent()
21562167
q = self.base_ring().order()
@@ -2160,8 +2171,8 @@ cdef class Polynomial(CommutativePolynomial):
21602171
return self
21612172

21622173
# We expect to succeed with greater than 1/2 probability,
2163-
# so if we try 1000 times and fail, there's a bug somewhere.
2164-
for _ in range(1000):
2174+
# so if we try 100 times and fail, there's a bug somewhere.
2175+
for _ in range(100):
21652176
# Sample a polynomial "uniformly" from R
21662177
# TODO: once #37118 has been merged, this can be made cleaner,
21672178
# as we will actually have access to uniform sampling.
@@ -2174,7 +2185,7 @@ cdef class Polynomial(CommutativePolynomial):
21742185

21752186
# Need to handle odd and even characteristic separately
21762187
if q % 2:
2177-
h = self.gcd(pow(T, (q-1)//2, self) - 1)
2188+
h = self.gcd(pow(T, (q**degree-1)//2, self) - 1)
21782189
else:
21792190
# Compute the trace of T with field of order 2^k
21802191
# sum T^(2^i) for i in range (degree * k)
@@ -2200,20 +2211,23 @@ cdef class Polynomial(CommutativePolynomial):
22002211
# If you are reaching this error, chances are there's a bug in the code.
22012212
raise AssertionError(f"no splitting of degree {degree} found for {self}")
22022213

2203-
def _any_irreducible_factor_squarefree(self, degree=None):
2214+
def _any_irreducible_factor_squarefree(self, degree=None, ext_degree=None):
22042215
"""
2205-
Helper function for any_irreducible_factor which computes
2216+
Helper function for :meth:`any_irreducible_factor` which computes
22062217
an irreducible factor from self, assuming the input is
22072218
squarefree.
22082219
22092220
Does this by first computing the distinct degree factorisations
22102221
of self and then finds a factor with Cantor-Zassenhaus
22112222
splitting.
22122223
2213-
If degree is not None, then only irreducible factors of degree
2214-
`degree` are searched for, otherwise the smallest degree factor
2224+
If ``degree`` is not ``None``, then only irreducible factors of degree
2225+
``degree`` are searched for, otherwise the smallest degree factor
22152226
is found.
22162227
2228+
If ``ext_degree`` is not ``None``, then only irreducible factors whose
2229+
degree divides ``ext_degree`` are returned.
2230+
22172231
EXAMPLES::
22182232
22192233
sage: # needs sage.rings.finite_rings
@@ -2262,10 +2276,15 @@ cdef class Polynomial(CommutativePolynomial):
22622276

22632277
# Otherwise we use the smallest possible d value
22642278
for (poly, d) in self._distinct_degree_factorisation_squarefree():
2265-
return poly._cantor_zassenhaus_split_to_irreducible(d)
2266-
raise ValueError(f"no irreducible factor could be computed from {self}")
2279+
if ext_degree is None:
2280+
return poly._cantor_zassenhaus_split_to_irreducible(d)
2281+
elif ZZ(d).divides(ext_degree):
2282+
return poly._cantor_zassenhaus_split_to_irreducible(d)
2283+
if d > ext_degree:
2284+
raise ValueError(f"no irreducible factor of degree {degree} dividing {ext_degree} could be computed from {self}")
2285+
raise AssertionError(f"no irreducible factor could be computed from {self}")
22672286

2268-
def any_irreducible_factor(self, degree=None, assume_squarefree=False, assume_distinct_deg=False):
2287+
def any_irreducible_factor(self, degree=None, assume_squarefree=False, assume_equal_deg=False, ext_degree=None):
22692288
"""
22702289
Return an irreducible factor of this polynomial.
22712290
@@ -2281,11 +2300,15 @@ cdef class Polynomial(CommutativePolynomial):
22812300
Used for polynomials over finite fields. If ``True``,
22822301
this polynomial is assumed to be squarefree.
22832302
2284-
- ``assume_distinct_deg`` (boolean) -- (default: ``False``).
2303+
- ``assume_equal_deg`` (boolean) -- (default: ``False``).
22852304
Used for polynomials over finite fields. If ``True``,
22862305
this polynomial is assumed to be the product of irreducible
22872306
polynomials of degree equal to ``degree``.
22882307
2308+
- ``ext_degree`` -- positive integer or ``None`` (default);
2309+
used for polynomials over finite fields. If not ``None`` only returns
2310+
irreducible factors of ``self`` whose degree divides ``ext_degree``.
2311+
22892312
EXAMPLES::
22902313
22912314
sage: # needs sage.rings.finite_rings
@@ -2328,9 +2351,9 @@ cdef class Polynomial(CommutativePolynomial):
23282351
sage: F = GF(163)
23292352
sage: R.<x> = F[]
23302353
sage: h = (x + 57) * (x + 98) * (x + 117) * (x + 145)
2331-
sage: h.any_irreducible_factor(degree=1, assume_distinct_deg=True) # random
2354+
sage: h.any_irreducible_factor(degree=1, assume_equal_deg=True) # random
23322355
x + 98
2333-
sage: h.any_irreducible_factor(assume_distinct_deg=True)
2356+
sage: h.any_irreducible_factor(assume_equal_deg=True)
23342357
Traceback (most recent call last):
23352358
...
23362359
ValueError: degree must be known if distinct degree factorisation is assumed
@@ -2359,7 +2382,7 @@ cdef class Polynomial(CommutativePolynomial):
23592382
if degree < 1:
23602383
raise ValueError(f"{degree = } must be positive")
23612384

2362-
if assume_distinct_deg and degree is None:
2385+
if assume_equal_deg and degree is None:
23632386
raise ValueError("degree must be known if distinct degree factorisation is assumed")
23642387

23652388
# When not working over a finite field, do the simple thing of factoring.
@@ -2397,27 +2420,30 @@ cdef class Polynomial(CommutativePolynomial):
23972420

23982421
# If we know the polynomial is square-free, we can start here
23992422
if assume_squarefree:
2400-
if assume_distinct_deg:
2423+
if assume_equal_deg:
24012424
return self._cantor_zassenhaus_split_to_irreducible(degree)
2402-
return self._any_irreducible_factor_squarefree(degree)
2425+
return self._any_irreducible_factor_squarefree(degree, ext_degree)
24032426

24042427
# Otherwise we compute the squarefree decomposition and check each
24052428
# polynomial for a root. If no poly has a root, we raise an error.
24062429
SFD = self.squarefree_decomposition()
24072430
SFD.sort()
24082431
for poly, _ in SFD:
24092432
try:
2410-
return poly._any_irreducible_factor_squarefree(degree)
2433+
return poly._any_irreducible_factor_squarefree(degree, ext_degree)
24112434
except ValueError:
24122435
pass
24132436

24142437
# If degree has been set, there could just be no factor of the desired degree
24152438
if degree:
24162439
raise ValueError(f"polynomial {self} has no irreducible factor of degree {degree}")
2440+
# If ext_degree has been set, then there may be no irreducible factor of degree dividing ext_degree
2441+
if ext_degree:
2442+
raise ValueError(f"polynomial {self} has no irreducible factor of degree dividing {ext_degree}")
24172443
# But if any degree is allowed then there should certainly be a factor if self has degree > 0
24182444
raise AssertionError(f"no irreducible factor was computed for {self}. Bug.")
24192445

2420-
def any_root(self, ring=None, degree=None, assume_squarefree=False, assume_distinct_deg=False):
2446+
def any_root(self, ring=None, degree=None, assume_squarefree=False, assume_equal_deg=False):
24212447
"""
24222448
Return a root of this polynomial in the given ring.
24232449
@@ -2437,14 +2463,26 @@ cdef class Polynomial(CommutativePolynomial):
24372463
finite fields. If ``True``, this polynomial is assumed to be
24382464
squarefree.
24392465
2440-
- ``assume_distinct_deg`` (bool) -- Used for polynomials over
2466+
- ``assume_equal_deg`` (bool) -- Used for polynomials over
24412467
finite fields. If ``True``, all factors of this polynomial
2442-
are assumed to have degree ``degree``.
2468+
are assumed to have degree ``degree``. Note that ``degree``
2469+
must be set.
24432470
24442471
.. WARNING::
24452472
24462473
Negative degree input will be deprecated. Instead use
2447-
``assume_distinct_deg``.
2474+
``assume_equal_deg``.
2475+
2476+
.. NOTE::
2477+
2478+
For finite fields, ``any_root()`` is non-deterministic when
2479+
finding linear roots of a polynomial over the base ring.
2480+
However, if ``degree`` is greater than one, or ``ring`` is an
2481+
extension of the base ring, then the root computed is found
2482+
by attempting to return a root after factorisation. Roots found
2483+
in this way are deterministic. This may change in the future.
2484+
For all other rings or fields, roots are found by first
2485+
fully-factoring ``self`` and the output is deterministic.
24482486
24492487
EXAMPLES::
24502488
@@ -2554,38 +2592,76 @@ cdef class Polynomial(CommutativePolynomial):
25542592
# When not working over a finite field, do the simple thing of factoring for
25552593
# roots and picking the first root. If none are available, raise an error.
25562594
from sage.categories.finite_fields import FiniteFields
2557-
if not self.base_ring() in FiniteFields():
2558-
rs = self.roots(ring=ring, multiplicities=False)
2559-
if rs:
2560-
return rs[0]
2561-
raise ValueError(f"polynomial {self} has no roots")
2595+
if self.base_ring() not in FiniteFields():
2596+
if ring not in FiniteFields():
2597+
rs = self.roots(ring=ring, multiplicities=False)
2598+
if rs:
2599+
return rs[0]
2600+
raise ValueError(f"polynomial {self} has no roots")
2601+
2602+
# Ensure that a provided ring is appropriate for the function. From the
2603+
# above we know it is either None or a finite field. When it's a finite
2604+
# field we ensure there's a coercion from the base ring to ring.
2605+
if ring is not None:
2606+
if ring.coerce_map_from(self.base_ring()) is None:
2607+
raise ValueError(f"no coercion map can be computed from {self.base_ring()} to {ring}")
25622608

25632609
# When the degree is none, we only look for a linear factor
25642610
if degree is None:
2565-
# if a ring is given try and coerce the polynomial into this ring
2566-
if ring is not None:
2611+
# When ring is None, we attempt to find a linear factor of self
2612+
if ring is None:
25672613
try:
2568-
self = self.change_ring(ring)
2614+
f = self.any_irreducible_factor(degree=1, assume_squarefree=assume_squarefree)
25692615
except ValueError:
2570-
raise(f"cannot coerce polynomial {self} to the new ring: {ring}")
2616+
raise ValueError(f"no root of polynomial {self} can be computed")
2617+
return - f[0] / f[1]
25712618

2572-
# try and find a linear irreducible polynomial from f to compute a root
2619+
# When we have a ring, then we can find an irreducible factor of degree `d` providing
2620+
# that d divides the degree of the extension from the base ring to the given ring
2621+
allowed_extension_degree = ring.degree() // self.base_ring().degree()
25732622
try:
2574-
f = self.any_irreducible_factor(degree=1, assume_squarefree=assume_squarefree)
2623+
f = self.any_irreducible_factor(assume_squarefree=assume_squarefree, ext_degree=allowed_extension_degree)
25752624
except ValueError:
2576-
raise ValueError(f"no root of polynomial {self} can be computed")
2577-
2578-
return - f[0] / f[1]
2625+
raise ValueError(f"no root of polynomial {self} can be computed over the ring {ring}")
2626+
# When d != 1 we then find the smallest extension
2627+
# TODO: What we should do here is compute some minimal
2628+
# extension F_ext = self.base_ring().extension(d, names="a") and find a
2629+
# root here and then coerce this root into the parent ring. This means we
2630+
# would work with the smallest possible extension.
2631+
# However, if we have some element of GF(p^k) and we try and coerce this to
2632+
# some element GF(p^(k*n)) this can fail, even though mathematically it
2633+
# should be fine.
2634+
# TODO: Additionally, if the above was solved, it would be faster to extend the base
2635+
# ring with the irreducible factor however, if the base ring is an extension
2636+
# then the type of self.base_ring().extension(f) is a Univariate Quotient Polynomial Ring
2637+
# and not a finite field.
2638+
2639+
# When f has degree one we simply return the roots
2640+
# TODO: should we write something fast for degree two using
2641+
# the quadratic formula?
2642+
if f.degree().is_one():
2643+
root = - f[0] / f[1]
2644+
return ring(root)
2645+
2646+
# TODO: The proper thing to do here would be to call
2647+
# return f.change_ring(ring).any_root()
2648+
# but as we cannot work in the minimal extension (see above) working
2649+
# in the extension for f.change_ring(ring).any_root() is almost always
2650+
# much much much slower than using the call for roots() which uses
2651+
# C library bindings for all finite fields.
2652+
# Until the coercion system for finite fields works better,
2653+
# this will be the most performant
2654+
return f.roots(ring, multiplicities=False)[0]
25792655

25802656
# The old version of `any_root()` allowed degree < 0 to indicate that the input polynomial
25812657
# had a distinct degree factorisation, we pass this to any_irreducible_factor as a bool and
25822658
# ensure that the degree is positive.
25832659
degree = ZZ(degree)
25842660
if degree < 0:
25852661
from sage.misc.superseded import deprecation
2586-
deprecation(37170, "negative ``degree`` will be disallowed. Instead use the bool `assume_distinct_deg`.")
2662+
deprecation(37170, "negative ``degree`` will be disallowed. Instead use the bool `assume_equal_deg`.")
25872663
degree = -degree
2588-
assume_distinct_deg = True
2664+
assume_equal_deg = True
25892665

25902666
# If a certain degree is requested, then we find an irreducible factor of degree `degree`
25912667
# use this to compute a field extension and return the generator as root of this polynomial
@@ -2594,7 +2670,7 @@ cdef class Polynomial(CommutativePolynomial):
25942670
try:
25952671
f = self.any_irreducible_factor(degree=degree,
25962672
assume_squarefree=assume_squarefree,
2597-
assume_distinct_deg=assume_distinct_deg)
2673+
assume_equal_deg=assume_equal_deg)
25982674
except ValueError:
25992675
raise ValueError(f"no irreducible factor of degree {degree} can be computed from {self}")
26002676

@@ -2617,13 +2693,17 @@ cdef class Polynomial(CommutativePolynomial):
26172693
# FiniteField type if the base field is a non-prime field,
26182694
# so this slower option is chosen to ensure the root is
26192695
# over explicitly a FiniteField type.
2620-
ring = self.base_ring().extension(f.degree(), names="a")
2621-
2622-
# Now we look for a linear root of this irreducible polynomial of degree `degree`
2623-
# over the user supplied ring or the extension we just computed. If the user sent
2624-
# a bad ring here of course there may be no root found.
2625-
f = f.change_ring(ring)
2626-
return f.any_root()
2696+
ring = self.base_ring().extension(degree, names="a")
2697+
2698+
# TODO: The proper thing to do here would be to call
2699+
# return f.change_ring(ring).any_root()
2700+
# but as we cannot work in the minimal extension (see above) working
2701+
# in the extension for f.change_ring(ring).any_root() is almost always
2702+
# much much much slower than using the call for roots() which uses
2703+
# C library bindings for all finite fields.
2704+
# Until the coercion system for finite fields works better,
2705+
# this will be the most performant
2706+
return f.roots(ring, multiplicities=False)[0]
26272707

26282708
def __truediv__(left, right):
26292709
r"""

src/sage/rings/polynomial/polynomial_quotient_ring.py

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -2098,7 +2098,7 @@ def _isomorphic_ring(self):
20982098
base_image = self.base_ring().modulus().change_ring(isomorphic_ring).any_root()
20992099
base_to_isomorphic_ring = self.base_ring().hom([isomorphic_ring(base_image)])
21002100
modulus = self.modulus().map_coefficients(base_to_isomorphic_ring)
2101-
gen = modulus.any_root(assume_squarefree=True, degree=1, assume_distinct_deg=True)
2101+
gen = modulus.any_root(assume_squarefree=True, degree=1, assume_equal_deg=True)
21022102

21032103
homspace = Hom(self, isomorphic_ring)
21042104
to_isomorphic_ring = homspace.__make_element_class__(SetMorphism)(homspace,

src/sage/schemes/elliptic_curves/cm.py

Lines changed: 11 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -233,6 +233,13 @@ def is_HCP(f, check_monic_irreducible=True):
233233
sage: all(not is_HCP(hilbert_class_polynomial(D) + 1)
234234
....: for D in srange(-4,-100,-1) if D.is_discriminant())
235235
True
236+
237+
Ensure that :issue:`37471` is fixed::
238+
239+
sage: from sage.schemes.elliptic_curves.cm import is_HCP
240+
sage: set_random_seed(297388353221545796156853787333338705098)
241+
sage: is_HCP(hilbert_class_polynomial(-55))
242+
-55
236243
"""
237244
zero = ZZ(0)
238245
# optional check that input is monic and irreducible
@@ -264,14 +271,15 @@ def is_HCP(f, check_monic_irreducible=True):
264271
# Compute X^p-X mod fp
265272
z = fp.parent().gen()
266273
r = pow(z, p, fp) - z
267-
d = r.gcd(fp).degree() # number of roots mod p
274+
r = r.gcd(fp)
275+
d = r.degree() # number of roots mod p
268276
if d == 0:
269277
continue
270-
if not fp.is_squarefree():
278+
if not r.is_squarefree():
271279
continue
272280
if d < h and d not in h2list:
273281
return zero
274-
jp = fp.any_root(degree=1, assume_squarefree=True, assume_distinct_deg=True)
282+
jp = r.any_root(degree=1, assume_squarefree=True, assume_equal_deg=True)
275283
E = EllipticCurve(j=jp)
276284
if E.is_supersingular():
277285
continue

0 commit comments

Comments
 (0)