Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
Download
241782 views
1
"""
2
3
(c) Copyright 2009-2010 Salman Baig and Chris Hall
4
5
This file is part of ELLFF
6
7
ELLFF is free software: you can redistribute it and/or modify
8
it under the terms of the GNU General Public License as published by
9
the Free Software Foundation, either version 3 of the License, or
10
(at your option) any later version.
11
12
ELLFF is distributed in the hope that it will be useful,
13
but WITHOUT ANY WARRANTY; without even the implied warranty of
14
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15
GNU General Public License for more details.
16
17
You should have received a copy of the GNU General Public License
18
along with this program. If not, see <http://www.gnu.org/licenses/>.
19
20
"""
21
22
# Wrapper code for ellff library
23
24
include "stdsage.pxi"
25
include "cdefs.pxi"
26
27
from sage.libs.ntl.ntl_ZZ_decl cimport *
28
from sage.libs.ntl.ntl_ZZ cimport ntl_ZZ
29
from sage.libs.ntl.ntl_ZZ_p_decl cimport *
30
from sage.libs.ntl.ntl_ZZ_pX_decl cimport *
31
from sage.libs.ntl.ntl_ZZ_pE_decl cimport *
32
from sage.libs.ntl.ntl_ZZ_pEX_decl cimport *
33
from sage.libs.ntl.ntl_ZZ_pX cimport ntl_ZZ_pX
34
from sage.libs.ntl.ntl_ZZ_pE cimport ntl_ZZ_pE
35
from sage.libs.ntl.ntl_ZZ_pEX cimport ntl_ZZ_pEX
36
from sage.libs.ntl.ntl_lzz_pX_decl cimport *
37
from sage.libs.ntl.ntl_lzz_pX cimport *
38
from sage.libs.ntl.ntl_lzz_pContext cimport *
39
from sage.libs.ntl.ntl_lzz_pContext import ntl_zz_pContext
40
from sage.rings.all import PolynomialRing, GF, ZZ
41
from sage.rings.ring import is_Field
42
from sage.rings.integer import is_Integer
43
from sage.rings.arith import is_prime
44
from sage.rings.integer cimport Integer
45
#from sage.rings.fraction_field_element import FractionFieldElement, is_FractionFieldElement
46
from sage.rings.polynomial.polynomial_element import is_Polynomial
47
from sage.schemes.elliptic_curves.ell_generic import is_EllipticCurve
48
from sage.schemes.elliptic_curves.constructor import EllipticCurve
49
from sage.structure.sage_object import SageObject
50
51
####################
52
# PSAGE imports
53
from psage.function_fields import is_FunctionField, FunctionField
54
from psage.function_fields.function_field import RationalFunctionField
55
from psage.function_fields.function_field_element import FunctionFieldElement_rational
56
####################
57
58
# import psage.ellff.euler_database as edb
59
60
####################
61
62
cdef extern from "ntl_wrap.h":
63
void ZZ_to_mpz(mpz_t *, ZZ_c *)
64
65
####################
66
67
cdef extern from "NTL/lzz_pE.h":
68
ctypedef struct zz_pEContext_c "struct zz_pEContext":
69
void (*restore)()
70
71
zz_pEContext_c* zz_pEContext_new "new zz_pEContext"(zz_pX_c p)
72
void zz_pEContext_delete "delete "(zz_pEContext_c *mem)
73
74
void zz_pEContext_restore(zz_pEContext_c *c)
75
76
ctypedef struct zz_pE_c "struct zz_pE":
77
pass
78
79
zz_pE_c* zz_pE_new "new zz_pE"()
80
void zz_pE_delete "delete "(zz_pE_c *mem)
81
void zz_pE_conv "conv"(zz_pE_c x, zz_pX_c a)
82
zz_pX_c zz_pE_rep "rep"(zz_pE_c z)
83
84
####################
85
86
cdef extern from "NTL/lzz_pEX.h":
87
ctypedef struct zz_pEX_c "struct zz_pEX":
88
pass
89
90
zz_pEX_c* zz_pEX_new "new zz_pEX"()
91
void zz_pEX_delete "delete "(zz_pEX_c *mem)
92
93
void zz_pEX_SetCoeff "SetCoeff"(zz_pEX_c x, long i, zz_pE_c a)
94
95
####################
96
97
cdef extern from "lzz_pEratX.h":
98
ctypedef struct zz_pEratX_c "zz_pEratX":
99
pass
100
101
zz_pEratX_c* zz_pEratX_new "new zz_pEratX"(zz_pEX_c num, zz_pEX_c den)
102
void zz_pEratX_delete "delete "(zz_pEratX_c *mem)
103
104
####################
105
106
cdef extern from "euler.h":
107
void __euler_table "euler_table"(long *table, long min_tau, long max_tau, int euler_deg) except+
108
void __twist_table "twist_table"(zz_pEX_c f, long *untwisted_table, long *twisted_table, long min_tau, long max_tau)
109
void __pullback_table "pullback_table"(zz_pEX_c finite_disc, zz_pEX_c infinite_disc, zz_pEratX_c f, long *base_table, long *pullback_table, long min_tau, long max_tau)
110
void __sum_table "sum_table"(ZZ_c b_n, long *table, long min_tau, long max_tau)
111
void __compute_b_n "compute_b_n"(ZZ_c b_n)
112
void __compute_c_n "compute_c_n"(ZZ_c **b, ZZ_c **c, int n)
113
114
####################
115
116
cdef extern from "jacobi.h":
117
void __jacobi_sum "jacobi_sum"(ZZ_c ***sum, long d)
118
119
####################
120
121
cdef extern from "ell_surface.h":
122
ctypedef struct ell_surfaceInfoT_c "ell_surfaceInfoT":
123
int sign
124
int deg_L
125
int constant_f
126
127
void ell_surface_init "ell_surface::init"(zz_pEratX_c a1, zz_pEratX_c a2, zz_pEratX_c a3, zz_pEratX_c a4, zz_pEratX_c a6)
128
129
ctypedef struct ell_surface_c "ell_surface":
130
int sign
131
int deg_L
132
int constant_f
133
134
void ell_surface_init "ell_surface::init"(zz_pEratX_c a1, zz_pEratX_c a2, zz_pEratX_c a3, zz_pEratX_c a4, zz_pEratX_c a6)
135
ell_surfaceInfoT_c* ell_surface_getSurfaceInfo "ell_surface::getSurfaceInfo"()
136
137
void ell_get_j_invariant "get_j_invariant"(ZZ_pEX_c j_num, ZZ_pEX_c j_den)
138
void ell_get_an "get_an"(ZZ_pEX_c a4, ZZ_pEX_c a6)
139
void ell_get_reduction "get_reduction"(ZZ_pEX_c **divisors, int finite)
140
void ell_get_disc "get_disc"(ZZ_pEX_c disc, int finite_f)
141
142
####################
143
144
cdef extern from "helper.h":
145
cdef void __get_modulus "get_modulus"(zz_pX_c pi_1, zz_pX_c pi_2, zz_pX_c a, int p, int d1, int d2)
146
cdef void init_NTL_ff(int p, int d, int precompute_inverses, int precompute_square_roots, int precompute_legendre_char, int precompute_pth_frobenius_map)
147
148
####################
149
150
cdef class _ellff_EllipticCurve_c:
151
cdef ZZ_c **b, **b_star, **c, **c_star, **__b, **__c
152
cdef int bc_len
153
154
cdef long **trace_tables, *trace_table_sizes
155
cdef int num_trace_tables
156
157
def __init__(self):
158
pass
159
160
def __cinit__(self):
161
# allocate memory for tables of traces. initially each
162
# table is empty (and will need to be allocated).
163
self.num_trace_tables = 15
164
self.trace_tables = \
165
<long**>sage_malloc(sizeof(long*)*self.num_trace_tables)
166
self.trace_table_sizes = \
167
<long*>sage_malloc(sizeof(long)*self.num_trace_tables)
168
for i in range(self.num_trace_tables):
169
self.trace_tables[i] = NULL
170
self.trace_table_sizes[i] = -1
171
172
# allocate memory for b and c coefficients (for L-function)
173
self.bc_len = 2*self.num_trace_tables
174
self.b = <ZZ_c **>sage_malloc(sizeof(ZZ_c*)*(self.bc_len+1))
175
self.c = <ZZ_c **>sage_malloc(sizeof(ZZ_c*)*(self.bc_len+1))
176
self.b_star = \
177
<ZZ_c **>sage_malloc(sizeof(ZZ_c*)*(self.bc_len+1))
178
self.c_star = \
179
<ZZ_c **>sage_malloc(sizeof(ZZ_c*)*(self.bc_len+1))
180
181
for n in range(self.bc_len+1):
182
self.b[n] = ZZ_new()
183
self.c[n] = ZZ_new()
184
self.b_star[n] = ZZ_new()
185
self.c_star[n] = ZZ_new()
186
187
ZZ_set_from_int(self.c[0], 1)
188
ZZ_set_from_int(self.c_star[0], 1)
189
190
def __dealloc__(self):
191
for i in range(self.num_trace_tables):
192
if self.trace_table_sizes[i] != -1:
193
sage_free(self.trace_tables[i])
194
sage_free(self.trace_tables)
195
sage_free(self.trace_table_sizes)
196
197
for i in range(self.bc_len+1):
198
ZZ_delete(self.b[i])
199
ZZ_delete(self.c[i])
200
ZZ_delete(self.b_star[i])
201
ZZ_delete(self.c_star[i])
202
203
sage_free(self.b)
204
sage_free(self.c)
205
sage_free(self.b_star)
206
sage_free(self.c_star)
207
208
def _build_euler_table(self, n):
209
self._allocate_table(n, self.q**n+1)
210
211
self._surface_init(n)
212
__euler_table(self.trace_tables[n], 0, self.q**n, self.d)
213
214
cdef Integer b_n = Integer(0)
215
for i in range(self.q**n+1):
216
b_n += self.trace_tables[n][i]
217
218
b_n._to_ZZ(self.b[n])
219
220
def _twist_euler_table(self, f, _ellff_EllipticCurve_c E, n):
221
cdef zz_pEX_c f_ = _to_zz_pEX_(f, self.phi, self.p)
222
223
assert self.trace_table_sizes[n] == self.q**n+1
224
assert E.trace_table_sizes[n] == self.q**n+1
225
__twist_table(f_, self.trace_tables[n], E.trace_tables[n],
226
0, self.q**n)
227
228
cdef Integer b_n = Integer(0)
229
for i in range(self.q**n+1):
230
b_n += E.trace_tables[n][i]
231
232
b_n._to_ZZ(E.b[n])
233
234
def _pullback_euler_table(self, f, _ellff_EllipticCurve_c E, n):
235
cdef zz_pEratX_c *f_
236
f_ = _to_EratX_(f, self.phi, self.p)
237
238
cdef zz_pEX_c fin_disc, inf_disc
239
fin_disc = _to_zz_pEX_(self.__finite_disc, self.phi, self.p)
240
inf_disc = _to_zz_pEX_(self.__infinite_disc, self.phi, self.p)
241
242
__pullback_table(fin_disc, inf_disc, f_[0],
243
self.trace_tables[n], E.trace_tables[n], 0, self.q**n)
244
245
cdef Integer b_n = Integer(0)
246
for i in range(self.q**n+1):
247
b_n += E.trace_tables[n][i]
248
249
b_n._to_ZZ(E.b[n])
250
ZZ_sub(E.b_star[n][0], E.b[n][0], self.b[n][0])
251
252
def _allocate_table(self, n, size):
253
if self.trace_table_sizes[n] == -1:
254
self.trace_tables[n] \
255
= <long*>sage_malloc(sizeof(long)*size)
256
self.trace_table_sizes[n] = size
257
elif self.trace_table_sizes[n] != size:
258
raise RuntimeError("table size is wrong")
259
260
def _deallocate_table(self, i):
261
assert self.trace_table_sizes[i] > 0
262
sage_free(self.trace_tables[i])
263
self.trace_table_sizes[i] = -1
264
265
def _num_trace_tables(self):
266
return self.num_trace_tables
267
268
def _trace_table_sizes(self, n):
269
return self.trace_table_sizes[n]
270
271
def _compute_c_n(self, n):
272
__compute_c_n(self.b, self.c, n)
273
274
def _compute_rel_c_n(self, n):
275
__compute_c_n(self.b_star, self.c_star, n)
276
277
def __common_functional_equation(self, euler_deg, deg_L, sign, optimal_deg):
278
cdef ZZ_c **b, **c
279
b = self.__b
280
c = self.__c
281
282
cdef ZZ_c fec_c, m_c
283
if euler_deg >= optimal_deg:
284
if deg_L % 2 == 1: m = self.q * sign
285
else: m = 1 * sign
286
287
assert not(deg_L % 2 == 0 and sign==1) or euler_deg*2 >= deg_L
288
289
ZZ_set_from_int(&m_c, m)
290
for n in range((deg_L+deg_L%2)/2, deg_L+1):
291
if n <= euler_deg:
292
if deg_L-n < n:
293
# compare what we calculated against f.e.
294
ZZ_mul(fec_c, m_c, c[deg_L-n][0])
295
if not ZZ_equal(c[n][0], fec_c):
296
e = 'Coefficient', n, \
297
'does not match. By hand:', \
298
ZZ_to_int(&c[n][0]), ', FE:', \
299
ZZ_to_int(&fec_c)
300
raise RuntimeError(e)
301
302
elif deg_L % 2 == 0 and n == deg_L/2:
303
assert sign == -1
304
# f.e. implies middle coefficient vanishes
305
ZZ_set_from_int(c[n], 0)
306
307
else:
308
ZZ_mul(c[n][0], m_c, c[deg_L-n][0])
309
310
ZZ_mul_long(m_c, m_c, self.q * self.q)
311
312
def _apply_functional_equation(self, euler_deg):
313
self.__b = self.b
314
self.__c = self.c
315
self.__common_functional_equation(euler_deg=euler_deg,
316
deg_L=self.deg_L,
317
sign=self.sign,
318
optimal_deg=self.optimal_deg)
319
self.b = self.__b
320
self.c = self.__c
321
322
# perform a further consistency check, if possible
323
if euler_deg > self.deg_L:
324
for n in range(self.deg_L+1,euler_deg+1):
325
if not ZZ_IsZero(self.c[n][0]):
326
return RuntimeError("fails a consistency check")
327
self.c_len = euler_deg
328
elif euler_deg >= self.optimal_deg:
329
self.c_len = self.deg_L
330
else:
331
self.c_len = euler_deg
332
333
def _apply_rel_functional_equation(self, euler_deg):
334
self.__b = self.b_star
335
self.__c = self.c_star
336
self.__common_functional_equation(euler_deg=euler_deg,
337
deg_L=self.rel_deg_L,
338
sign=self.rel_sign,
339
optimal_deg=self.rel_optimal_deg)
340
self.b_star = self.__b
341
self.c_star = self.__c
342
343
# perform a further consistency check, if possible
344
if euler_deg > self.rel_deg_L:
345
for n in range(self.rel_deg_L+1,euler_deg+1):
346
if not ZZ_IsZero(self.c[n][0]):
347
return RuntimeError("fails a consistency check")
348
self.c_star_len = euler_deg
349
elif euler_deg >= self.rel_optimal_deg:
350
self.c_star_len = self.rel_deg_L
351
else:
352
self.c_star_len = euler_deg
353
354
def _build_L_function(self, int euler_deg):
355
if euler_deg < self.optimal_deg: deg_L = euler_deg
356
else: deg_L = self.deg_L
357
358
cdef Integer tmp
359
v = []
360
for n in range(deg_L+1):
361
tmp = PY_NEW(Integer)
362
ZZ_to_mpz(&tmp.value, self.c[n])
363
v.append(tmp)
364
365
R = PolynomialRing(ZZ, 'T')
366
return R(v)
367
368
def _build_relative_L_function(self, int euler_deg):
369
if euler_deg < self.rel_optimal_deg: rel_deg_L = euler_deg
370
else: rel_deg_L = self.rel_deg_L
371
372
cdef Integer tmp
373
v = []
374
for n in range(rel_deg_L+1):
375
tmp = PY_NEW(Integer)
376
ZZ_to_mpz(&tmp.value, self.c_star[n])
377
v.append(tmp)
378
379
R = PolynomialRing(ZZ, 'T')
380
return R(v)
381
382
def __euler_table(self, n):
383
q = self.q
384
if self.trace_table_sizes[n] == -1:
385
raise RuntimeError("table is empty")
386
assert self.trace_table_sizes[n] == q**n+1
387
388
v = []
389
for i in range(0, q**n+1):
390
v.append(self.trace_tables[n][i])
391
return v
392
393
def __set_euler_table(self, n, table, force=False):
394
q = self.q
395
if self.trace_table_sizes[n] != -1:
396
if not force:
397
raise RuntimeError("run with force=True to force an overwrite")
398
self._deallocate_table(n)
399
400
self._allocate_table(n, q**n+1)
401
for i in range(q**n+1):
402
self.trace_tables[n][i] = table[i]
403
404
cdef Integer b_n = Integer(0)
405
for i in range(self.q**n+1):
406
b_n += self.trace_tables[n][i]
407
408
b_n._to_ZZ(self.b[n])
409
410
def _b(self, n):
411
cdef Integer tmp
412
tmp = PY_NEW(Integer)
413
ZZ_to_mpz(&tmp.value, self.b[n])
414
return tmp
415
416
def _c(self, n):
417
cdef Integer tmp
418
tmp = PY_NEW(Integer)
419
ZZ_to_mpz(&tmp.value, self.c[n])
420
return tmp
421
422
def _b_star(self, n):
423
cdef Integer tmp
424
tmp = PY_NEW(Integer)
425
ZZ_to_mpz(&tmp.value, self.b_star[n])
426
return tmp
427
428
def _c_star(self, n):
429
cdef Integer tmp
430
tmp = PY_NEW(Integer)
431
ZZ_to_mpz(&tmp.value, self.c_star[n])
432
return tmp
433
434
435
class ellff_EllipticCurve(_ellff_EllipticCurve_c,SageObject):
436
r"""
437
The \class{ellff_EllipticCurve} class represents an ELLFF
438
elliptic curve.
439
"""
440
441
def __init__(self, field, ainvs):
442
r"""
443
Create the ellff elliptic curve with invariants
444
\code{ainvs}, which is a list of 5 \emph{field elements}
445
$a_1$, $a_2$, $a_3$, $a_4$, and $a_6$.
446
447
INPUT:
448
- field -- ground field
449
- ainvs -- a list of 5 integers
450
451
EXAMPLES::
452
453
sage: import psage
454
sage: K.<t> = psage.FunctionField(GF(5))
455
sage: E = psage.ellff_EllipticCurve(K,[0,0,0,t^2,t+1]); E
456
<class 'psage.ellff.ellff.ellff_EllipticCurve'>
457
"""
458
459
_ellff_EllipticCurve_c.__init__(self)
460
461
if not is_Field(field):
462
raise TypeError, "field must be a field."
463
464
if not is_FunctionField(field):
465
raise TypeError, "field must be of type FunctionField"
466
467
if not isinstance(field, RationalFunctionField):
468
raise NotImplementedError, "K must be the function field of P^1."
469
470
# TODO: make this into a check for K=F_q(t)
471
# true/false: the following words *only if* K =F_q(C)
472
self.K = field
473
self.R = self.K.maximal_order()
474
475
F = self.F = self.K.constant_field()
476
p = self.p = self.F.characteristic()
477
d = self.d = self.F.degree()
478
q = self.q = p**d
479
480
self._E = EllipticCurve(field, ainvs)
481
self.L_function_calculated = False;
482
self.relative_L_function_calculated = False;
483
484
if not isinstance(ainvs, list) and len(ainvs) == 5:
485
raise TypeError, "ainvs must be a list of length 5."
486
487
cdef int sign, deg_L, constant_f
488
sign, deg_L, constant_f, phi = self._surface_init(1)
489
490
self.sign = sign
491
self.deg_L = deg_L
492
self.constant_f = constant_f
493
self.phi = phi
494
self.__calc_optimal_deg()
495
self.euler_deg = 0
496
self.rel_euler_deg = 0
497
498
def _retrieve_invariants(self):
499
r"""
500
Convert the invariants associated to self (local reduction
501
information for the finite and infinite models, the
502
$a$ invariants (only $a_4$ and $a_6$ for now), the $j$-invariant, and
503
the finite and infinite discriminant) from elements of Sage's
504
constant field to elements of ELLFF's constant field via a
505
field embedding.
506
507
EXAMPLES::
508
509
sage: import psage
510
sage: K.<t> = psage.FunctionField(GF(11))
511
sage: E = psage.ellff_EllipticCurve(K,[0,0,0,(t+1)*(t+2),t^2+1])
512
sage: E._retrieve_invariants()
513
514
"""
515
516
phi, F2 = _ellff_field_embedding(self.F, self.p, self.d, 1)
517
518
# respective finite fields for sage and ellff
519
sage_F = self.F
520
ellff_F = F2
521
522
cdef ZZ_pEX_c **divisors
523
divisors = <ZZ_pEX_c**>sage_malloc(sizeof(ZZ_pEX_c*)*9)
524
for i in range(9):
525
divisors[i] = ZZ_pEX_new()
526
527
# get information about bad reduction away from infinity
528
ell_get_reduction(divisors, 1)
529
self._finite_reduction = []
530
for i in range(9):
531
self._finite_reduction.append(
532
self.R(__from_ZZ_pEX(sage_F, ellff_F, divisors[i][0], self.R)))
533
534
# self._finite_reduction[0] = M_sp
535
# self._finite_reduction[1] = M_ns
536
# self._finite_reduction[2] = I^*
537
# self._finite_reduction[3] = II
538
# self._finite_reduction[4] = II^*
539
# self._finite_reduction[5] = III
540
# self._finite_reduction[6] = III^*
541
# self._finite_reduction[7] = IV
542
# self._finite_reduction[8] = IV^*
543
544
# convert into course information about reduction
545
self._finite_M = 1
546
self._finite_A = 1
547
for i in range(0,2):
548
self._finite_M *= self._finite_reduction[i]
549
for i in range(2,9):
550
self._finite_A *= self._finite_reduction[i]
551
552
# get information about reduction for minimal about infinity
553
ell_get_reduction(divisors, 0)
554
self._infinite_reduction = []
555
for i in range(9):
556
self._infinite_reduction.append(
557
self.R(__from_ZZ_pEX(sage_F, ellff_F, divisors[i][0], self.R)))
558
559
# self._infinite_reduction[0] = M_sp
560
# self._infinite_reduction[1] = M_ns
561
# self._infinite_reduction[2] = I^*
562
# self._infinite_reduction[3] = II
563
# self._infinite_reduction[4] = II^*
564
# self._infinite_reduction[5] = III
565
# self._infinite_reduction[6] = III^*
566
# self._infinite_reduction[7] = IV
567
# self._infinite_reduction[8] = IV^*
568
569
# convert into course information about reduction
570
self._infinite_M = 1
571
self._infinite_A = 1
572
for i in range(0,2):
573
self._infinite_M *= self._infinite_reduction[i]
574
for i in range(2,9):
575
self._infinite_A *= self._infinite_reduction[i]
576
577
cdef ZZ_pEX_c *a4 = ZZ_pEX_new(), *a6 = ZZ_pEX_new()
578
579
ell_get_an(a4[0], a6[0])
580
self.a4 = self.R(__from_ZZ_pEX(self.F, F2, a4[0], self.R))
581
self.a6 = self.R(__from_ZZ_pEX(self.F, F2, a6[0], self.R))
582
583
cdef ZZ_pEX_c *j_num = ZZ_pEX_new(), *j_den = ZZ_pEX_new()
584
585
ell_get_j_invariant(j_num[0], j_den[0])
586
self.j = self.R(__from_ZZ_pEX(self.F, F2, j_num[0], self.R)) / self.R(__from_ZZ_pEX(self.F, F2, j_den[0], self.R))
587
588
cdef ZZ_pEX_c *finite_disc = ZZ_pEX_new()
589
cdef ZZ_pEX_c *infinite_disc = ZZ_pEX_new()
590
ell_get_disc(finite_disc[0], 1)
591
ell_get_disc(infinite_disc[0], 0)
592
593
self.__finite_disc = self.R(__from_ZZ_pEX(self.F, F2, finite_disc[0], self.R))
594
self.__infinite_disc = self.R(__from_ZZ_pEX(self.F, F2, infinite_disc[0], self.R))
595
596
# clean up allocated memory
597
for i in range(9):
598
ZZ_pEX_delete(divisors[i])
599
sage_free(divisors)
600
601
ZZ_pEX_delete(a4)
602
ZZ_pEX_delete(a6)
603
ZZ_pEX_delete(j_num)
604
ZZ_pEX_delete(j_den)
605
ZZ_pEX_delete(finite_disc)
606
ZZ_pEX_delete(infinite_disc)
607
608
def _surface_init(self, n):
609
r"""
610
Initialize the elliptic surface associated to self over $F_{q^n}$.
611
612
INPUT:
613
614
- n -- an integer
615
616
OUTPUT:
617
618
- info.sign -- the sign of the functional equation
619
- info.deg_L -- the degree of the L-function
620
- info.constant_f -- True (1) if curve is constant and has no places of bad reduction
621
- phi -- an embedding of Sage's F_q to ELLFF's F_q
622
623
EXAMPLES::
624
625
sage: import psage
626
sage: K.<t> = psage.FunctionField(GF(11))
627
sage: E = psage.ellff_EllipticCurve(K,[0,0,0,(t+1)*(t+2),t^2+1])
628
sage: E._surface_init(2)
629
[1, 4, 0, Ring morphism:
630
From: Finite Field of size 11
631
To: Finite Field in c of size 11^2
632
Defn: 1 |--> 1]
633
634
"""
635
636
p = self.p
637
d = self.d
638
639
phi, F2 = _ellff_field_embedding(self.F, p, d, n)
640
641
# convert coefficients to format acceptible for ellff library
642
cdef zz_pEratX_c *a1_, *a2_, *a3_, *a4_, *a6_
643
644
a1_ = _to_EratX_(self._E.a1(), phi, p)
645
a2_ = _to_EratX_(self._E.a2(), phi, p)
646
a3_ = _to_EratX_(self._E.a3(), phi, p)
647
a4_ = _to_EratX_(self._E.a4(), phi, p)
648
a6_ = _to_EratX_(self._E.a6(), phi, p)
649
650
ell_surface_init(a1_[0], a2_[0], a3_[0], a4_[0], a6_[0])
651
652
zz_pEratX_delete(a1_)
653
zz_pEratX_delete(a2_)
654
zz_pEratX_delete(a3_)
655
zz_pEratX_delete(a4_)
656
zz_pEratX_delete(a6_)
657
658
init_NTL_ff(p, d*n, 0, 0, 0, 0)
659
660
cdef ell_surfaceInfoT_c *info
661
info = ell_surface_getSurfaceInfo()
662
663
# TODO: add context save (and later restore)
664
# cdef ell_surfaceContext_c *context
665
666
if n == 1:
667
self._retrieve_invariants()
668
669
return [info.sign, info.deg_L, info.constant_f, phi]
670
671
def __calc_optimal_deg(self):
672
r"""
673
Calculate and assign the optimal degree to compute the
674
$L$-function of self, i.e. what is the maximum degree of the
675
extension of $\mathbb{F}_q$ needed to determine the
676
$L$-function.
677
678
EXAMPLES::
679
680
sage: import psage
681
sage: K.<t> = psage.FunctionField(GF(11))
682
sage: E = psage.ellff_EllipticCurve(K,[0,0,0,(t+1)*(t+2),t^2+1])
683
sage: E.__calc_optimal_deg()
684
685
"""
686
687
if self.deg_L % 2 == 1:
688
self.optimal_deg = (self.deg_L-1)/2;
689
elif self.sign == 1:
690
self.optimal_deg = self.deg_L/2
691
else:
692
self.optimal_deg = self.deg_L/2 - 1
693
694
def quadratic_twist(self, f, tables=False, euler_deg=0, force=False, verbose=False):
695
r"""
696
Create a new elliptic curve which is a quadratic twist of this
697
curve. If specified, calculates tables of Euler factors for
698
twisted curve using tables of this curve.
699
700
INPUT:
701
f -- twisting polynomial (twist by extension K(sqrt{f})/K)
702
tables -- if True, computes Euler factors of twist.
703
(default: False)
704
euler_deg -- if positive and tables==True, a bound on the
705
degree of the Euler factors to compute. if zero, becomes
706
the optimal Euler degree for computing the L-function.
707
(default: 0)
708
force -- if True, make tables of Euler factors for this
709
curve, if necessary, as well as for new curve.
710
otherwise, use brute force.
711
(default: False)
712
verbose -- if True, outputs progress information.
713
(default: False)
714
715
EXAMPLES::
716
717
sage: import psage.ellff.ellff as ellff
718
sage: K.<t> = FunctionField(GF(5))
719
sage: E = ellff.ellff_EllipticCurve(K,[0,0,0,t^2,t+1]); E
720
<class 'ellff.ellff_EllipticCurve'>
721
sage: E_twist = E.quadratic_twist(t^2+1); E_twist
722
<class 'ellff.ellff_EllipticCurve'>
723
sage: E_twist.a4
724
t^8 + 2*t^6 + t^4
725
sage: E_twist.a6
726
t^9 + t^8 + 3*t^7 + 3*t^6 + 3*t^5 + 3*t^4 + t^3 + t^2
727
728
"""
729
730
q = self.q
731
732
self._surface_init(1)
733
734
E = ellff_EllipticCurve(self.K, [0, 0, 0, self.a4*f*f, self.a6*f*f*f])
735
736
cdef ntl_ZZ b_n = PY_NEW(ntl_ZZ)
737
if tables:
738
if euler_deg == 0:
739
euler_deg = E.optimal_deg
740
741
for n in range(1,euler_deg+1):
742
if force and self._trace_table_sizes(n) == -1:
743
if verbose:
744
print "rebuilding own Euler table (n = ", n, ")"
745
self._build_euler_table(n)
746
747
if self._trace_table_sizes(n) != -1:
748
if verbose:
749
print "twisting old table into new (n = ", n, ")"
750
E._surface_init(n)
751
E._allocate_table(n, q**n+1)
752
self._twist_euler_table(f, E, n)
753
754
else:
755
if verbose:
756
print "building new Euler table from scratch (n = ", n, ")"
757
E._build_euler_table(n)
758
759
return E
760
761
def pullback(self, f, tables=False, euler_deg=-1, verbose=False):
762
r"""
763
Create a new elliptic curve which is a pullback of this
764
curve. If specified, calculates tables of Euler factors for
765
pullback curve using tables of this curve. Will also
766
compute 'relative L-function', that is, the quotient
767
L(T,E_pullback)/L(T,E).
768
769
INPUT:
770
f -- Non-constant rational function. Gives rise to field
771
embedding F_q(f)-->F_q(t). Treat current curve as defined
772
over F_q(f) and pullback along field extension F_q(t)/F_q(t)
773
774
tables -- if True, computes Euler factors of pullback. When
775
base curve has tables of Euler factors, also computes
776
'relative b_n' (whose 'relative c_n' are coefficients of
777
the 'relative L-function').
778
(default: False)
779
780
euler_deg -- if non-negative and tables==True, equals
781
bound on the degree of Euler factors to compute. if negative,
782
becomes the optimal degree.
783
(default: None)
784
785
verbose -- if True, outputs progress information.
786
(default: False)
787
788
EXAMPLES::
789
790
sage: import psage.ellff.ellff as ellff
791
sage: K.<t> = FunctionField(GF(5))
792
sage: E = ellff.ellff_EllipticCurve(K,[0,0,0,t^2,t+1]); E
793
<class 'ellff.ellff_EllipticCurve'>
794
sage: E_pullback = E.pullback(t^3); E_pullback
795
<class 'ellff.ellff_EllipticCurve'>
796
sage: E_pullback.a4
797
t^4
798
sage: E_pullback.a6
799
t^3 + t^2
800
sage: E_pullback = E.pullback(t^3, tables=True, verbose=True); E_pullback
801
rebuilding own Euler table (n = 1 )
802
pulling back old table (n = 1 )
803
<class 'ellff.ellff_EllipticCurve'>
804
805
"""
806
807
q = self.q
808
809
self._surface_init(1)
810
811
new_a4 = self.a4.numerator().subs(f) / self.a4.denominator().subs(f)
812
new_a6 = self.a6.numerator().subs(f) / self.a6.denominator().subs(f)
813
E = ellff_EllipticCurve(self.K, [0, 0, 0, new_a4, new_a6])
814
815
E.rel_deg_L = E.deg_L - self.deg_L
816
E.rel_sign = E.sign / self.sign
817
E.rel_optimal_deg = E.optimal_deg - self.optimal_deg
818
819
cdef ntl_ZZ b_n = PY_NEW(ntl_ZZ)
820
if tables:
821
if euler_deg < 0:
822
euler_deg = E.rel_optimal_deg
823
E.rel_euler_deg = 0
824
825
for n in range(1,euler_deg+1):
826
if self._trace_table_sizes(n) == -1:
827
if verbose:
828
print "rebuilding own Euler table (n = ", n, ")"
829
self._build_euler_table(n)
830
831
if self._trace_table_sizes(n) != -1:
832
if verbose:
833
print "pulling back old table (n = ", n, ")"
834
E._surface_init(n)
835
E._allocate_table(n, q**n+1)
836
self._pullback_euler_table(f, E, n)
837
838
if E.rel_euler_deg == n-1:
839
E.rel_euler_deg = n
840
841
else:
842
if verbose:
843
print "building new Euler table (n = ", n, ")"
844
E._build_euler_table(n)
845
846
return E
847
848
def relative_L_function(self, verbose=False):
849
r"""
850
If not done so already, calculate the relative L-function
851
(or part of it) of the elliptic curve.
852
853
INPUT:
854
euler_deg -- if non-negative, a bound on the degree of the
855
Euler factors to compute. if negative, becomes the optimal
856
Euler degree for computing the L-function.
857
(default: -1)
858
859
EXAMPLES::
860
861
sage: import psage.ellff.ellff as ellff
862
sage: K.<t> = FunctionField(GF(5))
863
sage: E = ellff.ellff_EllipticCurve(K,[0,0,0,t^2,t+1]); E
864
<class 'ellff.ellff_EllipticCurve'>
865
sage: E_pullback = E.pullback(t^3, tables=True, verbose=True); E_pullback
866
rebuilding own Euler table (n = 1 )
867
pulling back old table (n = 1 )
868
<class 'ellff.ellff_EllipticCurve'>
869
sage: E_pullback.relative_L_function()
870
125*T^3 - 5*T^2 - T + 1
871
872
"""
873
874
if self.relative_L_function_calculated == True:
875
return self._relative_L_function
876
877
# retrieve useful constants
878
p = self.p
879
d = self.d
880
q = self.q = p**d
881
882
rel_deg_L = self.rel_deg_L
883
rel_sign = self.rel_sign
884
rel_optimal_deg = self.rel_optimal_deg
885
rel_euler_deg = self.rel_euler_deg
886
887
# determine relation between rel_euler_deg and rel_optimal_deg
888
if rel_optimal_deg > rel_euler_deg:
889
print "rel_euler_deg (= ", rel_euler_deg, " ) ", \
890
"is too small to compute entire L-function; ", \
891
"need rel_euler_deg >= ", rel_optimal_deg
892
893
if rel_euler_deg >= self._num_trace_tables():
894
raise RuntimeError("rel_euler_deg is too large")
895
896
# calculate necessary tables
897
for n in range(1, rel_euler_deg+1):
898
self._compute_rel_c_n(n)
899
900
# take functional equation into account, if possible
901
self._apply_rel_functional_equation(self.rel_euler_deg)
902
903
# produce output polynomial
904
self._relative_L_function = self._build_relative_L_function(rel_euler_deg)
905
906
if rel_euler_deg >= self.rel_optimal_deg:
907
self.relative_L_function_calculated = True
908
909
return self._relative_L_function
910
911
def L_function(self, int euler_deg=-1, force=False, verbose=False):
912
r"""
913
If not done so already, calculate the L-function (or part of it)
914
of the elliptic curve.
915
916
INPUT:
917
euler_deg -- if non-negative, a bound on the degree of the
918
Euler factors to compute. if negative, becomes the optimal
919
Euler degree for computing the L-function.
920
(default: 0)
921
force -- if True, forces recalculation of each table of
922
Euler factors used. if False, tables are only
923
constructed as needed.
924
(default: False)
925
926
EXAMPLES::
927
928
sage: import psage.ellff.ellff as ellff
929
sage: K.<t> = FunctionField(GF(5))
930
sage: E = ellff.ellff_EllipticCurve(K,[0,0,0,t^2,t+1]); E
931
<class 'ellff.ellff_EllipticCurve'>
932
sage: E.L_function()
933
625*T^4 + 125*T^3 + 5*T + 1
934
935
sage: E_pullback = E.pullback(t^3, tables=True, verbose=True); E_pullback
936
rebuilding own Euler table (n = 1 )
937
pulling back old table (n = 1 )
938
<class 'ellff.ellff_EllipticCurve'>
939
sage: E_pullback.L_function()
940
78125*T^7 + 1
941
942
sage: E_twist = E.quadratic_twist(t^2+1,tables=True,force=True); E_twist
943
<class 'ellff.ellff_EllipticCurve'>
944
sage: E_twist.L_function()
945
1953125*T^9 - 312500*T^8 - 31250*T^7 + 14750*T^6 - 10*T^5 - 2*T^4 + 118*T^3 - 10*T^2 - 4*T + 1
946
947
"""
948
949
if self.L_function_calculated == True:
950
return self._L_function
951
952
if self.constant_f == True:
953
raise ValueError, "Elliptic curve must be non-constant or have at least one place of bad reduction"
954
955
# retrieve useful constants
956
p = self.p
957
d = self.d
958
q = self.q = p**d
959
960
deg_L = self.deg_L
961
sign = self.sign
962
optimal_deg = self.optimal_deg
963
964
# determine relation between euler_deg and optimal_deg
965
if euler_deg < 0:
966
euler_deg = optimal_deg
967
elif deg_L > euler_deg*2:
968
print "euler_deg (= ", euler_deg, " ) ", \
969
"is too small to compute entire L-function; ", \
970
"need euler_deg >= ", optimal_deg
971
972
if euler_deg >= self._num_trace_tables():
973
raise RuntimeError("euler_deg is too large")
974
975
# calculate necessary tables
976
for n in range(1, euler_deg+1):
977
if force and self._trace_table_sizes(n) != -1:
978
self._deallocate_table(n)
979
980
if self._trace_table_sizes(n) == -1:
981
if verbose:
982
print "building euler table (n = ", n, ")"
983
self._allocate_table(n, q**n+1)
984
self._build_euler_table(n)
985
elif verbose:
986
print "already have euler table (n = ", n, ")"
987
988
self._compute_c_n(n)
989
990
if force or euler_deg >= self.euler_deg:
991
self.euler_deg = euler_deg
992
993
# take functional equation into account, if possible
994
self._apply_functional_equation(euler_deg)
995
996
# produce output polynomial
997
self._L_function = self._build_L_function(euler_deg)
998
999
if euler_deg >= optimal_deg:
1000
self.L_function_calculated = True
1001
1002
return self._L_function
1003
1004
def _euler_table(self, n):
1005
r"""
1006
Returns a python list representing the table of Euler factors
1007
1008
1009
INPUT:
1010
n -- which degree
1011
1012
EXAMPLES::
1013
1014
sage: import psage.ellff.ellff as ellff
1015
sage: K.<t> = FunctionField(GF(5))
1016
sage: E = ellff.ellff_EllipticCurve(K,[0,0,0,t^2,t+1]); E
1017
<class 'ellff.ellff_EllipticCurve'>
1018
sage: E._euler_table(1)
1019
Traceback (most recent call last):
1020
...
1021
RuntimeError: table is empty
1022
sage: E_pullback = E.pullback(t^3, tables=True, verbose=True); E_pullback
1023
rebuilding own Euler table (n = 1 )
1024
pulling back old table (n = 1 )
1025
<class 'ellff.ellff_EllipticCurve'>
1026
sage: E._euler_table(1)
1027
[0, 2, 3, -2, 2, 0]
1028
sage: E_pullback._euler_table(1)
1029
[0, 2, 1, -1, 2, 0]
1030
sage: E_twist = E.quadratic_twist(t^2+1, tables=True, force=True, verbose=True); E_twist
1031
twisting old table into new (n = 1 )
1032
twisting old table into new (n = 2 )
1033
twisting old table into new (n = 3 )
1034
twisting old table into new (n = 4 )
1035
<class 'ellff.ellff_EllipticCurve'>
1036
sage: E._euler_table(1)
1037
[0, 2, 3, -2, 2, 0]
1038
sage: E._euler_table(2)
1039
[-10, -6, -1, -6, -6, -4, 1, 4, -1, 1, -4, -1, 6, 1, -1, -4, -1, 6, 1, -1, -4, 1, 4, -1, 1, 0]
1040
sage: E._euler_table(3)[0:10]
1041
[0, -22, -18, 22, -22, 18, -12, -2, 3, 7]
1042
1043
sage: E_twist._euler_table(1)
1044
[0, -2, 0, 0, -2, 0]
1045
sage: E_twist._euler_table(2)
1046
[0, -6, 0, 0, -6, -4, -1, 4, -1, -1, -4, 1, -6, -1, 1, -4, 1, -6, -1, 1, -4, -1, 4, -1, -1, 0]
1047
sage: E_twist._euler_table(3)[0:10]
1048
[0, 22, 0, 0, 22, 18, 12, 2, 3, 7]
1049
1050
"""
1051
return self.__euler_table(n)
1052
1053
def _set_euler_table(self, n, table, force=False):
1054
r"""
1055
Sets the euler table of self over $F_{q^n}$ to table if table
1056
not already computed. If euler table has already been
1057
computed, it can be overwritten with force. Will also compute
1058
the corresponding exponential coefficients of the
1059
$L$-function, $b_n$.
1060
1061
INPUT:
1062
1063
- n -- the degree of the extenstion of F_q
1064
- table -- a table of euler factors
1065
- force -- a boolean that forces overwriting of existing euler table
1066
1067
EXAMPLES::
1068
1069
sage: import psage
1070
sage: K.<t> = psage.FunctionField(GF(11))
1071
sage: E = psage.ellff_EllipticCurve(K,[0,0,0,(t+1)*(t+2),t^2+1])
1072
sage: E._euler_table(1)
1073
Traceback (most recent call last):
1074
...
1075
RuntimeError: table is empty
1076
sage: E.L_function()
1077
14641*T^4 + 1
1078
sage: E._euler_table(1)
1079
[-4, -2, 1, -4, 3, -6, 3, 5, 4, 0, 0, 0]
1080
sage: et = E._euler_table(1)
1081
sage: E._set_euler_table(1,et)
1082
Traceback (most recent call last):
1083
...
1084
RuntimeError: run with force=True to force an overwrite
1085
sage: E._set_euler_table(1,table=et,force=True)
1086
1087
"""
1088
return self.__set_euler_table(n, table, force)
1089
1090
def _save_euler_table(self, n, verbose=False):
1091
r"""
1092
Save the euler table for self over the degree n extension of
1093
$\mathbb{F}_q$ to disk. This is currently implemented with
1094
sage.database.db, which uses ZODB. If self is the versal
1095
j-curve, it stores the table in the database
1096
1097
SAGE_ROOT/data/jcurve_euler_tables .
1098
1099
Otherwise, the tables are stored in the `user` table
1100
1101
SAGE_ROOT/data/local_euler_tables .
1102
1103
It currently doesn't check if the table already is stored; it
1104
merely writes over it in that case. This should eventually be
1105
implemented using MongoDB.
1106
1107
INPUT:
1108
1109
- n -- the degree of the extension of F_q
1110
1111
EXAMPLES::
1112
1113
sage: import psage
1114
sage: K.<t> = psage.FunctionField(GF(11))
1115
sage: E = psage.ellff_EllipticCurve(K,[0,0,0,-27*t/(t-1728),54*t/(t-1728)])
1116
sage: E._build_euler_table(1)
1117
sage: E._euler_table(1)
1118
[0, 0, 4, -6, 3, 5, 1, -2, 4, -2, 3, 1]
1119
sage: E._build_euler_table(2)
1120
sage: E._build_euler_table(3)
1121
sage: E._euler_table(1)
1122
[0, 0, 4, -6, 3, 5, 1, -2, 4, -2, 3, 1]
1123
sage: E._save_euler_table(1)
1124
sage: E._save_euler_table(2)
1125
sage: E._save_euler_table(3)
1126
1127
"""
1128
raise NotImplementedError( "ZODB was removed from Sage 5.7. Some other backend should be used here. See the port of Conway polynomials for a possible way to go" )
1129
# edb._save_euler_table(self, n, verbose)
1130
1131
def _load_euler_table(self, n, force=False, verbose=False):
1132
r"""
1133
Load the euler table for self over the degree n extension of
1134
$\mathbb{F}_q$ to disk. If self is the versal j-curve, the
1135
table is pulled from
1136
1137
SAGE_ROOT/data/jcurve_euler_tables .
1138
1139
Otherwise, the table is pulled from the `user` table
1140
1141
SAGE_ROOT/data/local_euler_tables .
1142
1143
This should eventually be implemented using MongoDB.
1144
1145
It currently doesn't check if the key exist. If the key
1146
doesn't exist, a RuntimeError is raised by
1147
sage.database.db. This RuntimeError should be sufficient, so
1148
key checking may not be necessary.
1149
1150
INPUT:
1151
1152
- n -- the degree of the extension of F_q
1153
1154
EXAMPLES::
1155
1156
sage: import psage
1157
sage: K.<t> = psage.FunctionField(GF(11))
1158
sage: E = psage.ellff_EllipticCurve(K,[0,0,0,-27*t/(t-1728),54*t/(t-1728)])
1159
sage: E._euler_table(1)
1160
Traceback (most recent call last):
1161
...
1162
RuntimeError: table is empty
1163
sage: E._load_euler_table(1)
1164
sage: E._euler_table(1)
1165
[0, 0, 4, -6, 3, 5, 1, -2, 4, -2, 3, 1]
1166
1167
"""
1168
raise NotImplementedError( "ZODB was removed from Sage 5.7. Some other backend should be used here. See the port of Conway polynomials for a possible way to go" )
1169
# edb._load_euler_table(self, n, force, verbose)
1170
1171
def M(self,fine=False):
1172
r"""
1173
Returns the divisor of multiplicative reduction.
1174
1175
INPUT:
1176
fine -- if True, returns pair of following outputs,
1177
one for split-multiplicative reduction and the
1178
other for non-split multiplicative reduction.
1179
1180
(default: False)
1181
1182
OUTPUT:
1183
( finite, infinite_degree )
1184
1185
finite -- polynomial in function field variable.
1186
the zeros are simple and precisely the places of
1187
(split or non-split) multiplicative reduction away
1188
from infinity.
1189
1190
infinite_degree -- 1 if there is (split or non-split)
1191
multiplicative reduction at infinity, and 0 otherwise.
1192
1193
EXAMPLES::
1194
1195
sage: import psage
1196
sage: K.<t> = psage.FunctionField(GF(11))
1197
sage: E = psage.ellff_EllipticCurve(K,[0,0,0,(t+1)*(t+2),t^2+1])
1198
sage: E.M()
1199
(t^6 + 9*t^5 + 4*t^4 + 8*t^3 + 8*t^2 + 3*t + 1, 0)
1200
1201
"""
1202
if fine == False:
1203
if self._infinite_M(0) == 0: return (self._finite_M, 1)
1204
else: return (self._finite_M, 0)
1205
1206
if self._infinite_reduction[0](0) == 0:
1207
p1 = (self._finite_reduction[0], 1)
1208
else:
1209
p1 = (self._finite_reduction[0], 0)
1210
1211
if self._infinite_reduction[1](0) == 0:
1212
p2 = (self._finite_reduction[1], 1)
1213
else:
1214
p2 = (self._finite_reduction[1], 0)
1215
1216
return (p1,p2)
1217
1218
def A(self,fine=False):
1219
r"""
1220
Returns the divisor additive reduction as a polynomial,
1221
for the finite part, and an integer.
1222
1223
INPUT:
1224
fine -- if True, returns tuple of following outputs,
1225
one for each type of bad reduction:
1226
I_n^*, II, II^*, III, III^*, IV, IV^*
1227
1228
(default: False)
1229
1230
OUTPUT:
1231
( finite, infinite_degree )
1232
1233
finite -- polynomial in function field variable.
1234
the zeros are simple and precisely the places of
1235
(special) additive reduction away from infinity.
1236
1237
infinite_degree -- 1 if there is (special) additive
1238
reduction at infinity, and 0 otherwise.
1239
1240
EXAMPLES::
1241
1242
sage: import psage
1243
sage: K.<t> = psage.FunctionField(GF(11))
1244
sage: E = psage.ellff_EllipticCurve(K,[0,0,0,(t+1)*(t+2),t^2+1])
1245
sage: E.A()
1246
(1, 1)
1247
1248
"""
1249
if fine == False:
1250
if self._infinite_A(0) == 0: return (self._finite_A, 1)
1251
else: return (self._finite_A, 0)
1252
1253
if self._infinite_reduction[2](0) == 0:
1254
p1 = (self._finite_reduction[2], 1)
1255
else:
1256
p1 = (self._finite_reduction[2], 0)
1257
1258
if self._infinite_reduction[3](0) == 0:
1259
p2 = (self._finite_reduction[3], 1)
1260
else:
1261
p2 = (self._finite_reduction[3], 0)
1262
1263
if self._infinite_reduction[4](0) == 0:
1264
p3 = (self._finite_reduction[4], 1)
1265
else:
1266
p3 = (self._finite_reduction[4], 0)
1267
1268
if self._infinite_reduction[5](0) == 0:
1269
p4 = (self._finite_reduction[5], 1)
1270
else:
1271
p4 = (self._finite_reduction[5], 0)
1272
1273
if self._infinite_reduction[6](0) == 0:
1274
p5 = (self._finite_reduction[6], 1)
1275
else:
1276
p5 = (self._finite_reduction[6], 0)
1277
1278
if self._infinite_reduction[7](0) == 0:
1279
p6 = (self._finite_reduction[7], 1)
1280
else:
1281
p6 = (self._finite_reduction[7], 0)
1282
1283
if self._infinite_reduction[8](0) == 0:
1284
p7 = (self._finite_reduction[8], 1)
1285
else:
1286
p7 = (self._finite_reduction[8], 0)
1287
1288
return (p1,p2,p3,p4,p5,p6,p7)
1289
1290
r"""
1291
def __repr__(self):
1292
a = self.ainvs()
1293
s = "y^2"
1294
if a[0] == -1:
1295
s += "- x*y "
1296
elif a[0] == 1:
1297
s += "+ x*y "
1298
elif a[0] != 0:
1299
s += "+ %s*x*y "%a[0]
1300
if a[2] == -1:
1301
s += " - y"
1302
elif a[2] == 1:
1303
s += "+ y"
1304
elif a[2] != 0:
1305
s += "+ %s*y"%a[2]
1306
s += " = x^3 "
1307
if a[1] == -1:
1308
s += "- x^2 "
1309
elif a[1] == 1:
1310
s += "+ x^2 "
1311
elif a[1] != 0:
1312
s += "+ %s*x^2 "%a[1]
1313
if a[3] == -1:
1314
s += "- x "
1315
elif a[3] == 1:
1316
s += "+ x "
1317
elif a[3] != 0:
1318
s += "+ %s*x "%a[3]
1319
if a[4] == -1:
1320
s += "-1"
1321
elif a[4] == 1:
1322
s += "+1"
1323
elif a[4] != 0:
1324
s += "+ %s"%a[4]
1325
s = s.replace("+ -","- ")
1326
return s
1327
"""
1328
1329
def conductor(self):
1330
r"""
1331
Return the conductor of this curve
1332
"""
1333
raise NotImplementedError
1334
1335
####################
1336
1337
## def my_func1(a1):
1338
## if is_FractionFieldElement(a1):
1339
## if is_Polynomial(a1.numerator()):
1340
## return 1
1341
## elif is_Polynomial(a1):
1342
## return 1
1343
## else:
1344
## raise TypeError, "Input must be rational function"
1345
1346
cdef zz_pE_c _to_zz_pE_(z, phi, p):
1347
cdef zz_pE_c* z_
1348
cdef ntl_zz_pX p_
1349
i_ = phi(z)
1350
if i_.category().object().degree() > 1:
1351
p_ = ntl_zz_pX(i_._vector_(), p)
1352
else:
1353
p_ = ntl_zz_pX([i_], p)
1354
z_ = zz_pE_new()
1355
zz_pE_conv(z_[0], p_.x)
1356
return z_[0]
1357
1358
## def _from_zz_pE(ntl_zz_pE z, F):
1359
## r"""
1360
## blah
1361
##
1362
## INPUT:
1363
## z -- element to embed.
1364
## F -- field to embed z into. must be isomorphic to field of z.
1365
##
1366
## (default: 0)
1367
##
1368
## EXAMPLES:
1369
## """
1370
## p = F.characteristic()
1371
## d = F.degree()
1372
##
1373
## # void zz_pE_conv "conv"(zz_pE_c x, zz_pX_c a)
1374
## # zz_pX_c zz_pE_rep "rep"(zz_pE_c z)
1375
1376
1377
cdef zz_pEX_c _to_zz_pEX_(f, phi, p):
1378
cdef zz_pEX_c* f_
1379
1380
f_ = zz_pEX_new()
1381
1382
# make sure f is really a polynomial
1383
if is_FunctionField(f.parent().fraction_field()):
1384
if f.denominator().degree() > 0:
1385
raise TypeError("f must be a polynomial")
1386
# it is, but it looks like a rational function
1387
f = f.numerator()
1388
1389
for i in range(f.degree()+1):
1390
zz_pEX_SetCoeff(f_[0], i, _to_zz_pE_(f.coeffs()[i], phi, p))
1391
1392
return f_[0]
1393
1394
1395
cdef zz_pEratX_c* _to_EratX_(r, phi, p):
1396
cdef zz_pEX_c n, d
1397
n = _to_zz_pEX_(r.numerator(), phi, p)
1398
d = _to_zz_pEX_(r.denominator(), phi, p)
1399
cdef zz_pEratX_c* r_
1400
# TODO: does this create a memory leak?
1401
r_ = zz_pEratX_new(n, d)
1402
return r_
1403
1404
1405
r"""
1406
blah
1407
1408
INPUT:
1409
blah --
1410
1411
(default: 0)
1412
1413
EXAMPLES:
1414
"""
1415
1416
cdef __from_ZZ_pE(sage_F, ellff_F, ZZ_pE_c z):
1417
phi = ellff_F.Hom(sage_F)([ellff_F.gen().minpoly().change_ring(sage_F).roots()[0][0]])
1418
1419
cdef ZZ_pX_c z_
1420
z_ = ZZ_pE_rep(z);
1421
d = ZZ_pX_deg(z_)
1422
v = []
1423
cdef long l
1424
for i in range(d+1):
1425
ZZ_conv_to_long(l, ZZ_p_rep(ZZ_pX_coeff(z_, i)))
1426
v.append(l)
1427
if (ZZ_pE_IsZero(z)):
1428
v = [0]
1429
if sage_F.degree() > 1:
1430
R = ZZ[ellff_F.variable_name()]
1431
v_ = R(v)
1432
return phi(ellff_F(v_))
1433
return sage_F(v[0])
1434
1435
cdef __from_ZZ_pEX(sage_F, ellff_F, ZZ_pEX_c f, R):
1436
r"""
1437
blah
1438
1439
INPUT:
1440
ring -- polynomial ring to embed in
1441
1442
(default: 0)
1443
1444
EXAMPLES:
1445
"""
1446
1447
d = ZZ_pEX_deg(f)
1448
cdef ZZ_pE_c c
1449
if not is_FunctionField(R.fraction_field()):
1450
raise TypeError("R.fraction_field() must be an order in a function field")
1451
# v = []
1452
# for i in range(d+1):
1453
# c = ZZ_pEX_coeff(f, i)
1454
# C = __from_ZZ_pE(sage_F, ellff_F, c)
1455
# v.append(C)
1456
#
1457
# return v
1458
s = R.gen()
1459
poly = R(0)
1460
for i in reversed(range(d+1)):
1461
c = ZZ_pEX_coeff(f, i)
1462
C = __from_ZZ_pE(sage_F, ellff_F, c)
1463
poly = s*poly + C
1464
return poly
1465
1466
def _ellff_field_embedding(F, p, d, n):
1467
r"""
1468
Return an embedding of F, an arbitrary finite field in Sage, to
1469
F2, a Sage finite field compatible with ELLFF. This embedding
1470
respects how Sage and ELLFF both treat relative extensions.
1471
1472
INPUT:
1473
- F -- a finite field in Sage
1474
- p -- the characteristic of F
1475
- d -- F_q = F_{p^d}
1476
- n -- F2 = F_{q^n}
1477
1478
OUTPUT:
1479
- phi -- the homomorphism F --> F2
1480
- F2 -- the ELLFF finite field as a Sage object
1481
1482
EXAMPLES::
1483
1484
sage: import psage
1485
sage: psage.ellff.ellff._ellff_field_embedding(GF(11), 11, 2, 3)
1486
[Ring morphism:
1487
From: Finite Field of size 11
1488
To: Finite Field in c of size 11^6
1489
Defn: 1 |--> 1, Finite Field in c of size 11^6]
1490
1491
"""
1492
cdef ntl_zz_pX pi_1 = ntl_zz_pX(modulus=p)
1493
cdef ntl_zz_pX pi_2 = ntl_zz_pX(modulus=p)
1494
cdef ntl_zz_pX alpha = ntl_zz_pX(modulus=p)
1495
1496
# get modulus from ellff library for representing F_q=F_{p^d},F_{q^n}
1497
__get_modulus(pi_1.x, pi_2.x, alpha.x, p, d, n)
1498
1499
# get embedding F_q-->F_{q^r}, r=relative_degree
1500
if d == 1:
1501
F2 = GF(p**pi_2.degree(), name='c', modulus=ZZ['x'](pi_2.list()))
1502
phi = F.Hom(F2)([F.gen().minpoly().change_ring(F2).roots()[0][0]])
1503
else:
1504
F1 = GF(p**pi_1.degree(), name='b', modulus=ZZ['x'](pi_1.list()))
1505
F2 = GF(p**pi_2.degree(), name='c', modulus=ZZ['x'](pi_2.list()))
1506
1507
phi1 = F.Hom(F1)([F.gen().minpoly().change_ring(F1).roots()[0][0]])
1508
phi2 = F1.Hom(F2)(ZZ['x'](alpha.list()))
1509
phi = phi1.post_compose(phi2)
1510
1511
return [phi, F2]
1512
1513
def jacobi_sum(p, n, d, verbose=False):
1514
r"""
1515
Returns a polynomial which represents an abstract Jacobi sum
1516
over the field F_q=F_{p^n}.
1517
1518
Suppose d|q-1 and let chi : F_q^* --> Z/d be the composition
1519
of x|-->x^{(q-1)/d} and a chosen isomorphism
1520
F_q^*/F_q^{*d} --> Z/d. let z1,z2 be independent generic
1521
elements satisfying z1^d=z2^d=1. The routine calculates the
1522
'abstract' Jacobi sum
1523
1524
J(d,q) = sum_{x!=0,1} z1^{chi(x)}*z2^{chi(1-x)}
1525
1526
and returns result as array A=(a_ij) where
1527
1528
a_ij = coefficient of z1^i*z2^j in J(d,q)
1529
1530
INPUT:
1531
1532
- d -- the abstract order of the characters; must divide q-1
1533
1534
OUTPUT:
1535
1536
A
1537
1538
EXAMPLES::
1539
1540
sage: import psage
1541
sage: psage.ellff.ellff.jacobi_sum(11,3,4)
1542
Traceback (most recent call last):
1543
...
1544
RuntimeError: d must divide q-1
1545
sage: psage.ellff.ellff.jacobi_sum(11,1,2)
1546
[[2, 2], [2, 3]]
1547
sage: psage.ellff.ellff.jacobi_sum(11,2,2)
1548
[[29, 30], [30, 30]]
1549
sage: psage.ellff.ellff.jacobi_sum(11,3,2)
1550
[[332, 332], [332, 333]]
1551
1552
"""
1553
1554
if (not is_Integer(p)) or (not is_prime(p)):
1555
raise RuntimeError("p must be prime")
1556
1557
if (not is_Integer(n)) or (n <= 0):
1558
raise RuntimeError("n must be a positive integer")
1559
1560
q = p ** n
1561
if (q-1) % d != 0:
1562
raise RuntimeError("d must divide q-1")
1563
1564
# calculate order of q mod d
1565
x = q % d
1566
f = 1
1567
while x != 1:
1568
x = (x*q) % d
1569
f = f+1
1570
if verbose:
1571
print "q = ", q, ", f = ", f
1572
1573
# allocate memory for Jacobi sum
1574
cdef ZZ_c ***sum
1575
sum = <ZZ_c***>sage_malloc(sizeof(ZZ_c**)*d)
1576
for i in range(d):
1577
sum[i] = <ZZ_c**>sage_malloc(sizeof(ZZ_c*)*d)
1578
for j in range(d):
1579
sum[i][j] = ZZ_new()
1580
1581
if verbose:
1582
print "setting up finite field"
1583
init_NTL_ff(p, n, 0, 0, 0, 0)
1584
1585
if verbose:
1586
print "calling library's jacobi_sum()"
1587
__jacobi_sum(sum, d)
1588
1589
# extract sum
1590
if verbose:
1591
print "extracting sum"
1592
cdef Integer tmp
1593
A = []
1594
for i in range(d):
1595
w = []
1596
for j in range(d):
1597
tmp = PY_NEW(Integer)
1598
ZZ_to_mpz(&tmp.value, sum[i][j])
1599
w.append(tmp)
1600
A.append(w)
1601
1602
# free memory we borrowed
1603
for i in range(d):
1604
for j in range(d):
1605
ZZ_delete(sum[i][j])
1606
sage_free(sum[i])
1607
sage_free(sum)
1608
1609
return A
1610
1611
def j_curve(K):
1612
r"""
1613
Returns the following elliptic curve over K:
1614
1615
y^2 = x^3 - 108*t/(t-1728)*x + 432*j/(j-1728).
1616
1617
INPUT:
1618
1619
- K -- function field F_q(t)
1620
1621
OUTPUT:
1622
1623
ellff_EllipticCurve
1624
1625
EXAMPLES::
1626
1627
sage: import psage
1628
sage: K.<t> = psage.FunctionField(GF(11))
1629
sage: E = psage.ellff.ellff.j_curve(K); E
1630
<class 'psage.ellff.ellff.ellff_EllipticCurve'>
1631
sage: E.a4, E.a6
1632
(2*t^4 + 5*t^3 + 6*t^2 + 9*t, 3*t^6 + 7*t^5 + 8*t^4 + 3*t^3 + 4*t^2 + 8*t)
1633
1634
"""
1635
1636
t = K.gens()[0]
1637
return ellff_EllipticCurve(K, [0,0,0,-108*t/(t-1728),432*t/(t-1728)])
1638
1639
1640