57
57
from sage .rings .power_series_ring import PowerSeriesRing
58
58
from . import hyperelliptic_generic
59
59
from sage .misc .cachefunc import cached_method
60
+ from sage .misc .superseded import deprecated_function_alias
60
61
from sage .matrix .constructor import identity_matrix , matrix
61
62
from sage .misc .functional import rank
62
63
from sage .misc .lazy_import import lazy_import
@@ -297,7 +298,7 @@ def frobenius_matrix(self, N=None, algorithm='hypellfrob'):
297
298
ValueError: In the current implementation, p must be greater than (2g+1)(2N-1) = 81
298
299
"""
299
300
if algorithm != 'hypellfrob' :
300
- raise ValueError ("Unknown algorithm" )
301
+ raise ValueError ("unknown algorithm" )
301
302
302
303
# By default, use precision enough to be able to compute the
303
304
# frobenius minimal polynomial
@@ -306,48 +307,41 @@ def frobenius_matrix(self, N=None, algorithm='hypellfrob'):
306
307
307
308
return self .frobenius_matrix_hypellfrob (N = N )
308
309
309
- def frobenius_polynomial_cardinalities (self , a = None ):
310
+ def _frobenius_polynomial_cardinalities (self , a = None ):
310
311
r"""
311
- Compute the charpoly of frobenius, as an element of `\ZZ[x]`,
312
- by computing the number of points on the curve over `g` extensions
313
- of the base field where `g` is the genus of the curve.
314
-
315
- .. WARNING::
316
-
317
- This is highly inefficient when the base field or the genus of the
318
- curve are large.
312
+ Helper method for :meth:`frobenius_polynomial`.
319
313
320
314
EXAMPLES::
321
315
322
316
sage: R.<t> = PolynomialRing(GF(37))
323
317
sage: H = HyperellipticCurve(t^5 + t + 2)
324
- sage: H.frobenius_polynomial_cardinalities( )
318
+ sage: H.frobenius_polynomial(algorithm='cardinalities' )
325
319
x^4 + x^3 - 52*x^2 + 37*x + 1369
326
320
327
321
A quadratic twist::
328
322
329
323
sage: H = HyperellipticCurve(2*t^5 + 2*t + 4)
330
- sage: H.frobenius_polynomial_cardinalities( )
324
+ sage: H.frobenius_polynomial(algorithm='cardinalities' )
331
325
x^4 - x^3 - 52*x^2 - 37*x + 1369
332
326
333
327
Curve over a non-prime field::
334
328
335
329
sage: K.<z> = GF(7**2)
336
330
sage: R.<t> = PolynomialRing(K)
337
331
sage: H = HyperellipticCurve(t^5 + z*t + z^2)
338
- sage: H.frobenius_polynomial_cardinalities( )
332
+ sage: H.frobenius_polynomial(algorithm='cardinalities' )
339
333
x^4 + 8*x^3 + 70*x^2 + 392*x + 2401
340
334
341
335
This method may actually be useful when ``hypellfrob`` does not work::
342
336
343
337
sage: K = GF(7)
344
338
sage: R.<t> = PolynomialRing(K)
345
339
sage: H = HyperellipticCurve(t^9 + t^3 + 1)
346
- sage: H.frobenius_polynomial_matrix (algorithm='hypellfrob ')
340
+ sage: H.frobenius_polynomial (algorithm='matrix ')
347
341
Traceback (most recent call last):
348
342
...
349
343
ValueError: In the current implementation, p must be greater than (2g+1)(2N-1) = 81
350
- sage: H.frobenius_polynomial_cardinalities( )
344
+ sage: H.frobenius_polynomial(algorithm='cardinalities' )
351
345
x^8 - 5*x^7 + 7*x^6 + 36*x^5 - 180*x^4 + 252*x^3 + 343*x^2 - 1715*x + 2401
352
346
"""
353
347
g = self .genus ()
@@ -372,37 +366,33 @@ def frobenius_polynomial_cardinalities(self, a=None):
372
366
373
367
return ZZ ['x' ](coeffs ).reverse ()
374
368
375
- def frobenius_polynomial_matrix (self , M = None , algorithm = 'hypellfrob' ):
369
+ def _frobenius_polynomial_matrix (self , M = None , algorithm = 'hypellfrob' ):
376
370
r"""
377
- Compute the charpoly of frobenius, as an element of `\ZZ[x]`,
378
- by computing the charpoly of the frobenius matrix.
379
-
380
- This is currently only supported when the base field is prime
381
- and large enough using the ``hypellfrob`` library.
371
+ Helper method for :meth:`frobenius_polynomial`.
382
372
383
373
EXAMPLES::
384
374
385
375
sage: R.<t> = PolynomialRing(GF(37))
386
376
sage: H = HyperellipticCurve(t^5 + t + 2)
387
- sage: H.frobenius_polynomial_matrix( )
377
+ sage: H.frobenius_polynomial(algorithm='matrix' )
388
378
x^4 + x^3 - 52*x^2 + 37*x + 1369
389
379
390
380
A quadratic twist::
391
381
392
382
sage: H = HyperellipticCurve(2*t^5 + 2*t + 4)
393
- sage: H.frobenius_polynomial_matrix( )
383
+ sage: H.frobenius_polynomial(algorithm='matrix' )
394
384
x^4 - x^3 - 52*x^2 - 37*x + 1369
395
385
396
386
Curves defined over larger prime fields::
397
387
398
388
sage: K = GF(49999)
399
389
sage: R.<t> = PolynomialRing(K)
400
390
sage: H = HyperellipticCurve(t^9 + t^5 + 1)
401
- sage: H.frobenius_polynomial_matrix( )
391
+ sage: H.frobenius_polynomial(algorithm='matrix' )
402
392
x^8 + 281*x^7 + 55939*x^6 + 14144175*x^5 + 3156455369*x^4 + 707194605825*x^3
403
393
+ 139841906155939*x^2 + 35122892542149719*x + 6249500014999800001
404
394
sage: H = HyperellipticCurve(t^15 + t^5 + 1)
405
- sage: H.frobenius_polynomial_matrix( ) # long time, 8s on a Corei7
395
+ sage: H.frobenius_polynomial(algorithm='matrix' ) # long time, 8s on a Corei7
406
396
x^14 - 76*x^13 + 220846*x^12 - 12984372*x^11 + 24374326657*x^10 - 1203243210304*x^9
407
397
+ 1770558798515792*x^8 - 74401511415210496*x^7 + 88526169366991084208*x^6
408
398
- 3007987702642212810304*x^5 + 3046608028331197124223343*x^4
@@ -440,52 +430,51 @@ def frobenius_polynomial_matrix(self, M=None, algorithm='hypellfrob'):
440
430
441
431
return ZZ ['x' ](f )
442
432
443
- def frobenius_polynomial_pari (self ):
433
+ def _frobenius_polynomial_pari (self ):
444
434
r"""
445
- Compute the charpoly of frobenius, as an element of `\ZZ[x]`,
446
- by calling the PARI function ``hyperellcharpoly``.
435
+ Helper method for :meth:`frobenius_polynomial`.
447
436
448
437
EXAMPLES::
449
438
450
439
sage: R.<t> = PolynomialRing(GF(37))
451
440
sage: H = HyperellipticCurve(t^5 + t + 2)
452
- sage: H.frobenius_polynomial_pari( )
441
+ sage: H.frobenius_polynomial(algorithm='pari' )
453
442
x^4 + x^3 - 52*x^2 + 37*x + 1369
454
443
455
444
A quadratic twist::
456
445
457
446
sage: H = HyperellipticCurve(2*t^5 + 2*t + 4)
458
- sage: H.frobenius_polynomial_pari( )
447
+ sage: H.frobenius_polynomial(algorithm='pari' )
459
448
x^4 - x^3 - 52*x^2 - 37*x + 1369
460
449
461
450
Slightly larger example::
462
451
463
452
sage: K = GF(2003)
464
453
sage: R.<t> = PolynomialRing(K)
465
454
sage: H = HyperellipticCurve(t^7 + 487*t^5 + 9*t + 1)
466
- sage: H.frobenius_polynomial_pari( )
455
+ sage: H.frobenius_polynomial(algorithm='pari' )
467
456
x^6 - 14*x^5 + 1512*x^4 - 66290*x^3 + 3028536*x^2 - 56168126*x + 8036054027
468
457
469
458
Curves defined over a non-prime field are supported as well::
470
459
471
460
sage: K.<a> = GF(7^2)
472
461
sage: R.<t> = PolynomialRing(K)
473
462
sage: H = HyperellipticCurve(t^5 + a*t + 1)
474
- sage: H.frobenius_polynomial_pari( )
463
+ sage: H.frobenius_polynomial(algorithm='pari' )
475
464
x^4 + 4*x^3 + 84*x^2 + 196*x + 2401
476
465
477
466
sage: K.<z> = GF(23**3)
478
467
sage: R.<t> = PolynomialRing(K)
479
468
sage: H = HyperellipticCurve(t^3 + z*t + 4)
480
- sage: H.frobenius_polynomial_pari( )
469
+ sage: H.frobenius_polynomial(algorithm='pari' )
481
470
x^2 - 15*x + 12167
482
471
483
472
Over prime fields of odd characteristic, `h` may be nonzero::
484
473
485
474
sage: K = GF(101)
486
475
sage: R.<t> = PolynomialRing(K)
487
476
sage: H = HyperellipticCurve(t^5 + 27*t + 3, t)
488
- sage: H.frobenius_polynomial_pari( )
477
+ sage: H.frobenius_polynomial(algorithm='pari' )
489
478
x^4 + 2*x^3 - 58*x^2 + 202*x + 10201
490
479
491
480
TESTS:
@@ -495,17 +484,54 @@ def frobenius_polynomial_pari(self):
495
484
sage: P.<x> = PolynomialRing(GF(3))
496
485
sage: u = x^10 + x^9 + x^8 + x
497
486
sage: C = HyperellipticCurve(u)
498
- sage: C.frobenius_polynomial_pari( )
487
+ sage: C.frobenius_polynomial(algorithm='pari' )
499
488
x^8 + 2*x^7 + 6*x^6 + 9*x^5 + 18*x^4 + 27*x^3 + 54*x^2 + 54*x + 81
500
489
"""
501
490
f , h = self .hyperelliptic_polynomials ()
502
491
return ZZ ['x' ](pari ([f , h ]).hyperellcharpoly ())
503
492
493
+ frobenius_polynomial_cardinalities = deprecated_function_alias (
494
+ 40528 , _frobenius_polynomial_cardinalities ,
495
+ replacement = 'frobenius_polynomial(algorithm="cardinalities")' ,
496
+ replacement_rst_doc = ':meth:`frobenius_polynomial(algorithm="cardinalities") <frobenius_polynomial>`' )
497
+ frobenius_polynomial_matrix = deprecated_function_alias (
498
+ 40528 , _frobenius_polynomial_matrix ,
499
+ replacement = 'frobenius_polynomial(algorithm="matrix")' ,
500
+ replacement_rst_doc = ':meth:`frobenius_polynomial(algorithm="matrix") <frobenius_polynomial>`' )
501
+ frobenius_polynomial_pari = deprecated_function_alias (
502
+ 40528 , _frobenius_polynomial_pari ,
503
+ replacement = 'frobenius_polynomial(algorithm="pari")' ,
504
+ replacement_rst_doc = ':meth:`frobenius_polynomial(algorithm="pari") <frobenius_polynomial>`' )
505
+
504
506
@cached_method
505
- def frobenius_polynomial (self ):
507
+ def frobenius_polynomial (self , algorithm = None , ** kwargs ):
506
508
r"""
507
509
Compute the charpoly of frobenius, as an element of `\ZZ[x]`.
508
510
511
+ INPUT:
512
+
513
+ - ``algorithm`` -- the algorithm to use, one of
514
+
515
+ - ``None`` (default) -- automatically select an algorithm.
516
+
517
+ - ``'pari'`` -- use the PARI function ``hyperellcharpoly``.
518
+
519
+ - ``'cardinalities'`` -- compute the number of points on the curve over `g`
520
+ extensions of the base field where `g` is the genus of the curve.
521
+
522
+ .. WARNING::
523
+
524
+ This is highly inefficient when the base field or the genus of the
525
+ curve are large.
526
+
527
+ - ``'matrix'`` -- compute the charpoly of the frobenius matrix.
528
+
529
+ This is currently only supported when the base field is prime
530
+ and large enough using the ``hypellfrob`` library.
531
+
532
+ - ``**kwargs`` -- additional keyword arguments passed to the
533
+ internal method. Undocumented feature, may be removed in the future.
534
+
509
535
EXAMPLES::
510
536
511
537
sage: R.<t> = PolynomialRing(GF(37))
@@ -578,6 +604,15 @@ def frobenius_polynomial(self):
578
604
sage: C.frobenius_polynomial()
579
605
x^8 + 2*x^7 + 6*x^6 + 9*x^5 + 18*x^4 + 27*x^3 + 54*x^2 + 54*x + 81
580
606
"""
607
+ if algorithm == "cardinalities" :
608
+ return self ._frobenius_polynomial_cardinalities (** kwargs )
609
+ elif algorithm == "matrix" :
610
+ return self ._frobenius_polynomial_matrix (** kwargs )
611
+ elif algorithm == "pari" :
612
+ return self ._frobenius_polynomial_pari (** kwargs )
613
+ elif algorithm is not None :
614
+ raise ValueError (f"unknown algorithm { algorithm } " )
615
+
581
616
K = self .base_ring ()
582
617
e = K .degree ()
583
618
q = K .cardinality ()
@@ -588,11 +623,11 @@ def frobenius_polynomial(self):
588
623
if (e == 1 and
589
624
q >= (2 * g + 1 )* (2 * self ._frobenius_coefficient_bound_charpoly ()- 1 ) and
590
625
h == 0 and f .degree () % 2 ):
591
- return self .frobenius_polynomial_matrix ( )
626
+ return self .frobenius_polynomial ( algorithm = 'matrix' )
592
627
elif q % 2 == 1 :
593
- return self .frobenius_polynomial_pari ( )
628
+ return self .frobenius_polynomial ( algorithm = 'pari' )
594
629
else :
595
- return self .frobenius_polynomial_cardinalities ( )
630
+ return self .frobenius_polynomial ( algorithm = 'cardinalities' )
596
631
597
632
def _points_fast_sqrt (self ):
598
633
"""
@@ -1031,7 +1066,7 @@ def count_points_exhaustive(self, n=1, naive=False):
1031
1066
a .append (self .cardinality_exhaustive (extension_degree = i ))
1032
1067
1033
1068
# let's not be too naive and compute the frobenius polynomial
1034
- f = self .frobenius_polynomial_cardinalities (a = a )
1069
+ f = self ._frobenius_polynomial_cardinalities (a = a )
1035
1070
return self .count_points_frobenius_polynomial (n = n , f = f )
1036
1071
1037
1072
def count_points_hypellfrob (self , n = 1 , N = None , algorithm = None ):
@@ -1116,7 +1151,7 @@ def count_points_hypellfrob(self, n=1, N=None, algorithm=None):
1116
1151
M = self .frobenius_matrix (N = N , algorithm = 'hypellfrob' )
1117
1152
return self .count_points_matrix_traces (n = n ,M = M ,N = N )
1118
1153
elif algorithm == 'charpoly' :
1119
- f = self .frobenius_polynomial_matrix (algorithm = 'hypellfrob' )
1154
+ f = self ._frobenius_polynomial_matrix (algorithm = 'hypellfrob' )
1120
1155
return self .count_points_frobenius_polynomial (n = n ,f = f )
1121
1156
else :
1122
1157
raise ValueError ("Unknown algorithm" )
0 commit comments