Line data Source code
1 : /* Copyright (C) 2010 The PARI group.
2 :
3 : This file is part of the PARI/GP package.
4 :
5 : PARI/GP is free software; you can redistribute it and/or modify it under the
6 : terms of the GNU General Public License as published by the Free Software
7 : Foundation; either version 2 of the License, or (at your option) any later
8 : version. It is distributed in the hope that it will be useful, but WITHOUT
9 : ANY WARRANTY WHATSOEVER.
10 :
11 : Check the License for details. You should have received a copy of it, along
12 : with the package; see the file 'COPYING'. If not, write to the Free Software
13 : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
14 :
15 : /********************************************************************/
16 : /** **/
17 : /** L functions of elliptic curves **/
18 : /** **/
19 : /********************************************************************/
20 : #include "pari.h"
21 : #include "paripriv.h"
22 :
23 : #define DEBUGLEVEL DEBUGLEVEL_ellanal
24 :
25 : struct baby_giant
26 : {
27 : GEN baby, giant, sum;
28 : GEN bnd, rbnd;
29 : };
30 :
31 : /* Generic Buhler-Gross algorithm */
32 :
33 : struct bg_data
34 : {
35 : GEN E, N; /* ell, conductor */
36 : GEN bnd; /* t_INT; will need all an for n <= bnd */
37 : ulong rootbnd; /* sqrt(bnd) */
38 : GEN an; /* t_VECSMALL: cache of an, n <= rootbnd */
39 : GEN p; /* t_VECSMALL: primes <= rootbnd */
40 : };
41 :
42 : typedef void bg_fun(void*el, GEN n, GEN a);
43 :
44 : /* a = a_n, where p = bg->pp[i] divides n, and lasta = a_{n/p}.
45 : * Call fun(E, N, a_N), for all N, n | N, P^+(N) <= p, a_N != 0,
46 : * i.e. assumes that fun accumulates a_N * w(N) */
47 :
48 : static void
49 1661632 : gen_BG_add(void *E, bg_fun *fun, struct bg_data *bg, GEN n, long i, GEN a, GEN lasta)
50 : {
51 1661632 : pari_sp av = avma;
52 : long j;
53 1661632 : ulong nn = itou_or_0(n);
54 1661632 : if (nn && nn <= bg->rootbnd) bg->an[nn] = itos(a);
55 :
56 1661632 : if (signe(a))
57 : {
58 418600 : fun(E, n, a);
59 418600 : j = 1;
60 : }
61 : else
62 1243032 : j = i;
63 3319701 : for(; j <= i; j++)
64 : {
65 2638230 : ulong p = bg->p[j];
66 2638230 : GEN nexta, pn = mului(p, n);
67 2638230 : if (cmpii(pn, bg->bnd) > 0) return;
68 1658069 : nexta = mulis(a, bg->an[p]);
69 1658069 : if (i == j && umodiu(bg->N, p)) nexta = subii(nexta, mului(p, lasta));
70 1658069 : gen_BG_add(E, fun, bg, pn, j, nexta, a);
71 1658069 : set_avma(av);
72 : }
73 : }
74 :
75 : static void
76 42 : gen_BG_init(struct bg_data *bg, GEN E, GEN N, GEN bnd)
77 : {
78 42 : bg->E = E;
79 42 : bg->N = N;
80 42 : bg->bnd = bnd;
81 42 : bg->rootbnd = itou(sqrtint(bnd));
82 42 : bg->p = primes_upto_zv(bg->rootbnd);
83 42 : bg->an = ellanQ_zv(E, bg->rootbnd);
84 42 : }
85 :
86 : static void
87 42 : gen_BG_rec(void *E, bg_fun *fun, struct bg_data *bg)
88 : {
89 42 : long i, j, lp = lg(bg->p)-1;
90 42 : GEN bndov2 = shifti(bg->bnd, -1);
91 42 : pari_sp av = avma, av2;
92 : GEN p;
93 : forprime_t S;
94 42 : (void)forprime_init(&S, utoipos(bg->p[lp]+1), bg->bnd);
95 42 : av2 = avma;
96 42 : if (DEBUGLEVEL)
97 0 : err_printf("1st stage, using recursion for p <= %ld\n", bg->p[lp]);
98 3605 : for (i = 1; i <= lp; i++)
99 : {
100 3563 : ulong pp = bg->p[i];
101 3563 : long ap = bg->an[pp];
102 3563 : gen_BG_add(E, fun, bg, utoipos(pp), i, stoi(ap), gen_1);
103 3563 : set_avma(av2);
104 : }
105 42 : if (DEBUGLEVEL) err_printf("2nd stage, looping for p <= %Ps\n", bndov2);
106 1157667 : while ( (p = forprime_next(&S)) )
107 : {
108 : long jmax;
109 1157667 : GEN ap = ellap(bg->E, p);
110 1157667 : pari_sp av3 = avma;
111 1157667 : if (!signe(ap)) continue;
112 :
113 578515 : jmax = itou( divii(bg->bnd, p) ); /* 2 <= jmax <= el->rootbound */
114 578515 : fun(E, p, ap);
115 9218090 : for (j = 2; j <= jmax; j++)
116 : {
117 8639575 : long aj = bg->an[j];
118 : GEN a, n;
119 8639575 : if (!aj) continue;
120 1352092 : a = mulis(ap, aj);
121 1352092 : n = muliu(p, j);
122 1352092 : fun(E, n, a);
123 1352092 : set_avma(av3);
124 : }
125 578515 : set_avma(av2);
126 578515 : if (abscmpii(p, bndov2) >= 0) break;
127 : }
128 42 : if (DEBUGLEVEL) err_printf("3nd stage, looping for p <= %Ps\n", bg->bnd);
129 1041467 : while ( (p = forprime_next(&S)) )
130 : {
131 1041425 : GEN ap = ellap(bg->E, p);
132 1041425 : if (!signe(ap)) continue;
133 520373 : fun(E, p, ap);
134 520373 : set_avma(av2);
135 : }
136 42 : set_avma(av);
137 42 : }
138 :
139 : /******************************************************************
140 : *
141 : * L functions of elliptic curves
142 : * Pascal Molin (molin.maths@gmail.com) 2014
143 : *
144 : ******************************************************************/
145 :
146 : struct lcritical
147 : {
148 : GEN h; /* real */
149 : long cprec; /* computation prec */
150 : long L; /* number of points */
151 : GEN K; /* length of series */
152 : long real;
153 : };
154 :
155 : static double
156 245 : logboundG0(long e, double aY)
157 : {
158 : double cla, loggam;
159 245 : cla = 1 + 1/sqrt(aY);
160 245 : if (e) cla = ( cla + 1/(2*aY) ) / (2*sqrt(aY));
161 245 : loggam = (e) ? M_LN2-aY : -aY + log( log( 1+1/aY) );
162 245 : return log(cla) + loggam;
163 : }
164 :
165 : static void
166 245 : param_points(GEN N, double Y, double tmax, long bprec, long *cprec, long *L,
167 : GEN *K, double *h)
168 : {
169 : double D, a, aY, X, logM;
170 245 : long d = 2, w = 1;
171 245 : tmax *= d;
172 245 : D = bprec * M_LN2 + M_PI/4*tmax + 2;
173 245 : *cprec = nbits2prec(ceil(D / M_LN2) + 5);
174 245 : a = 2 * M_PI / sqrt(gtodouble(N));
175 245 : aY = a * cos(M_PI/2*Y);
176 245 : logM = 2*M_LN2 + logboundG0(w+1, aY) + tmax * Y * M_PI/2;
177 245 : *h = ( 2 * M_PI * M_PI / 2 * Y ) / ( D + logM );
178 245 : X = log( D / a);
179 245 : *L = ceil( X / *h);
180 245 : *K = ceil_safe(dbltor( D / a ));
181 245 : }
182 :
183 : static GEN
184 245 : vecF2_lk(GEN E, GEN K, GEN rbnd, GEN Q, GEN sleh, long prec)
185 : {
186 : pari_sp av;
187 245 : long l, L = lg(K)-1;
188 245 : GEN a = ellanQ_zv(E, itos(gel(K,1)));
189 245 : GEN S = cgetg(L+1, t_VEC);
190 :
191 13776 : for (l = 1; l <= L; l++) gel(S,l) = cgetr(prec);
192 245 : av = avma;
193 13776 : for (l = 1; l <= L; l++)
194 : {
195 : GEN e1, Sl, z, zB;
196 13531 : long aB, b, A, B, Kl = itou(gel(K,l));
197 : pari_sp av2;
198 : /* FIXME: could reduce prec here (useful for large prec) */
199 13531 : e1 = gel(Q, l);
200 13531 : Sl = real_0(prec);;
201 : /* baby-step giant step */
202 13531 : B = A = rbnd[l];
203 13531 : z = powersr(e1, B); zB = gel(z, B+1);
204 13531 : av2 = avma;
205 314174 : for (aB = A*B; aB >= 0; aB -= B)
206 : {
207 300643 : GEN s = real_0(prec); /* could change also prec here */
208 38032519 : for (b = B; b > 0; --b)
209 : {
210 37731876 : long k = aB+b;
211 37731876 : if (k <= Kl && a[k]) s = addrr(s, mulsr(a[k], gel(z, b+1)));
212 37731876 : if (gc_needed(av2, 1)) (void)gc_all(av2, 2, &s, &Sl);
213 : }
214 300643 : Sl = addrr(mulrr(Sl, zB), s);
215 : }
216 13531 : affrr(mulrr(Sl, gel(sleh,l)), gel(S, l)); /* to avoid copying all S */
217 13531 : set_avma(av);
218 : }
219 245 : return S;
220 : }
221 :
222 : /* Return C, C[i][j] = Q[j]^i, i = 1..nb */
223 : static void
224 0 : baby_init(struct baby_giant *bb, GEN Q, GEN bnd, GEN rbnd, long prec)
225 : {
226 0 : long i, j, l = lg(Q);
227 : GEN R, C;
228 0 : C = cgetg(l,t_VEC);
229 0 : for (i = 1; i < l; ++i) gel(C, i) = powersr(gel(Q, i), rbnd[i]);
230 0 : R = cgetg(l,t_VEC);
231 0 : for (i = 1; i < l; ++i)
232 : {
233 0 : gel(R, i) = cgetg(rbnd[i]+1, t_VEC);
234 0 : gmael(R, i, 1) = rtor(gmael(C, i, 2), prec);
235 0 : for (j = 2; j <= rbnd[i]; j++) gmael(R, i, j) = stor(0, prec);
236 : }
237 0 : bb->baby = C; bb->giant = R;
238 0 : bb->bnd = bnd; bb->rbnd = rbnd;
239 0 : }
240 :
241 : static long
242 245 : baby_size(GEN rbnd, long Ks, long prec)
243 : {
244 245 : long i, s, m, l = lg(rbnd);
245 13776 : for (s = 0, i = 1; i < l; ++i)
246 13531 : s += rbnd[i];
247 245 : m = 2*s*prec + 3*l + s;
248 245 : if (DEBUGLEVEL)
249 0 : err_printf("ellL1: BG_add: %ld words, ellan: %ld words\n", m, Ks);
250 245 : return m;
251 : }
252 :
253 : static void
254 0 : ellL1_add(void *E, GEN n, GEN a)
255 : {
256 0 : pari_sp av = avma;
257 0 : struct baby_giant *bb = (struct baby_giant*) E;
258 0 : long j, l = lg(bb->giant);
259 0 : for (j = 1; j < l; j++)
260 0 : if (cmpii(n, gel(bb->bnd,j)) <= 0)
261 : {
262 0 : ulong r, q = uabsdiviu_rem(n, bb->rbnd[j], &r);
263 0 : GEN giant = gel(bb->giant, j), baby = gel(bb->baby, j);
264 0 : affrr(addrr(gel(giant, q+1), mulri(gel(baby, r+1), a)), gel(giant, q+1));
265 0 : set_avma(av);
266 0 : } else break;
267 0 : }
268 :
269 : static GEN
270 0 : vecF2_lk_bsgs(GEN E, GEN bnd, GEN rbnd, GEN Q, GEN sleh, GEN N, long prec)
271 : {
272 : struct bg_data bg;
273 : struct baby_giant bb;
274 0 : long k, L = lg(bnd)-1;
275 : GEN S;
276 0 : baby_init(&bb, Q, bnd, rbnd, prec);
277 0 : gen_BG_init(&bg, E, N, gel(bnd,1));
278 0 : gen_BG_rec((void*) &bb, ellL1_add, &bg);
279 0 : S = cgetg(L+1, t_VEC);
280 0 : for (k = 1; k <= L; ++k)
281 : {
282 0 : pari_sp av = avma;
283 0 : long j, g = rbnd[k];
284 0 : GEN giant = gmael(bb.baby, k, g+1), Sl = gmael(bb.giant, k, g);
285 0 : for (j = g-1; j >=1; j--) Sl = addrr(mulrr(Sl, giant), gmael(bb.giant,k,j));
286 0 : gel(S, k) = gc_uptoleaf(av, mulrr(gel(sleh,k), Sl));
287 : }
288 0 : return S;
289 : }
290 :
291 : static long
292 13531 : _sqrt(GEN x) { pari_sp av = avma; return gc_long(av, itou(sqrtint(x))); }
293 :
294 : static GEN
295 245 : vecF(struct lcritical *C, GEN E)
296 : {
297 245 : pari_sp av = avma;
298 245 : long prec = C->cprec, Ks = itos_or_0(C->K), L = C->L, l;
299 245 : GEN N = ellQ_get_N(E), PiN;
300 245 : GEN e = mpexp(C->h), elh = powersr(e, L-1), Q, bnd, rbnd, vec;
301 :
302 245 : PiN = divrr(Pi2n(1,prec), sqrtr_abs(itor(N, prec)));
303 245 : setsigne(PiN, -1); /* - 2Pi/sqrt(N) */
304 245 : bnd = gpowers0(invr(e), L-1, C->K); /* bnd[i] = K exp(-(i-1)h) */
305 245 : rbnd = cgetg(L+1, t_VECSMALL);
306 245 : Q = cgetg(L+1, t_VEC);
307 13776 : for (l = 1; l <= L; l++)
308 : {
309 13531 : gel(bnd,l) = ceil_safe(gel(bnd,l));
310 13531 : rbnd[l] = _sqrt(gel(bnd,l)) + 1;
311 13531 : gel(Q, l) = mpexp(mulrr(PiN, gel(elh, l)));
312 : }
313 245 : if (Ks && baby_size(rbnd, Ks, prec) > (Ks>>1))
314 245 : vec = vecF2_lk(E, bnd, rbnd, Q, elh, prec);
315 : else
316 0 : vec = vecF2_lk_bsgs(E, bnd, rbnd, Q, elh, N, prec);
317 245 : return gc_upto(av, vec);
318 : }
319 :
320 : /* Lambda function by Fourier inversion. vec is a grid, t a scalar or t_SER */
321 : static GEN
322 273 : glambda(GEN t, GEN vec, GEN h, long real, long prec)
323 : {
324 273 : GEN z, r, e = gexp(gmul(mkcomplex(gen_0,h), t), prec);
325 273 : long n = lg(vec)-1, i;
326 :
327 273 : r = real == 1? gmul2n(real_i(gel(vec, 1)), -1): gen_0;
328 273 : z = real == 1? e: gmul(powIs(3), e);
329 : /* FIXME: summing backward may be more stable */
330 16268 : for (i = 2; i <= n; i++)
331 : {
332 15995 : r = gadd(r, real_i(gmul(gel(vec,i), z)));
333 15995 : if (i < n) z = gmul(z, e);
334 : }
335 273 : return gmul(mulsr(4, h), r);
336 : }
337 :
338 : static GEN
339 245 : Lpoints(struct lcritical *C, GEN e, double tmax, long bprec)
340 : {
341 245 : double h = 0, Y = .97;
342 245 : GEN N = ellQ_get_N(e);
343 245 : param_points(N, Y, tmax, bprec, &C->cprec, &C->L, &C->K, &h);
344 245 : C->real = ellrootno_global(e);
345 245 : C->h = rtor(dbltor(h), C->cprec);
346 245 : return vecF(C, e);
347 : }
348 :
349 : static GEN
350 273 : Llambda(GEN vec, struct lcritical *C, GEN t, long prec)
351 : {
352 273 : GEN lambda = glambda(gprec_w(t, C->cprec), vec, C->h, C->real, C->cprec);
353 273 : return gprec_w(lambda, prec);
354 : }
355 :
356 : /* 2*(2*Pi)^(-s)*gamma(s)*N^(s/2); */
357 : static GEN
358 273 : ellgammafactor(GEN N, GEN s, long prec)
359 : {
360 273 : GEN c = gpow(divrr(gsqrt(N,prec), Pi2n(1,prec)), s, prec);
361 273 : return gmul(gmul2n(c,1), ggamma(s, prec));
362 : }
363 :
364 : static GEN
365 273 : ellL1_eval(GEN e, GEN vec, struct lcritical *C, GEN t, long prec)
366 : {
367 273 : GEN g = ellgammafactor(ellQ_get_N(e), gaddgs(gmul(gen_I(),t), 1), prec);
368 273 : return gdiv(Llambda(vec, C, t, prec), g);
369 : }
370 :
371 : static GEN
372 273 : ellL1_der(GEN e, GEN vec, struct lcritical *C, GEN t, long der, long prec)
373 : {
374 273 : GEN r = polcoef_i(ellL1_eval(e, vec, C, t, prec), der, 0);
375 273 : r = gmul(r,powIs(C->real == 1 ? -der: 1-der));
376 273 : return gmul(real_i(r), mpfact(der));
377 : }
378 :
379 : GEN
380 231 : ellL1(GEN E, long r, long bitprec)
381 : {
382 231 : pari_sp av = avma;
383 : struct lcritical C;
384 231 : long prec = nbits2prec(bitprec);
385 : GEN e, vec, t;
386 231 : if (r < 0)
387 7 : pari_err_DOMAIN("ellL1", "derivative order", "<", gen_0, stoi(r));
388 224 : e = ellanal_globalred(E, NULL);
389 224 : if (r == 0 && ellrootno_global(e) < 0) { set_avma(av); return gen_0; }
390 210 : vec = Lpoints(&C, e, 0., bitprec);
391 210 : t = r ? scalarser(gen_1, 0, r): zeroser(0, 0);
392 210 : setvalser(t, 1);
393 210 : return gc_upto(av, ellL1_der(e, vec, &C, t, r, prec));
394 : }
395 :
396 : GEN
397 35 : ellanalyticrank(GEN E, GEN eps, long bitprec)
398 : {
399 35 : pari_sp av = avma, av2;
400 35 : long prec = nbits2prec(bitprec);
401 : struct lcritical C;
402 : pari_timer ti;
403 : GEN e, vec;
404 : long rk;
405 35 : if (DEBUGLEVEL) timer_start(&ti);
406 35 : if (!eps)
407 35 : eps = real2n(-bitprec/2+1, DEFAULTPREC);
408 : else
409 0 : if (typ(eps) != t_REAL) {
410 0 : eps = gtofp(eps, DEFAULTPREC);
411 0 : if (typ(eps) != t_REAL) pari_err_TYPE("ellanalyticrank", eps);
412 : }
413 35 : e = ellanal_globalred(E, NULL);
414 35 : vec = Lpoints(&C, e, 0., bitprec);
415 35 : if (DEBUGLEVEL) timer_printf(&ti, "init L");
416 35 : av2 = avma;
417 35 : for (rk = C.real>0 ? 0: 1; ; rk += 2)
418 28 : {
419 : GEN Lrk;
420 63 : GEN t = rk ? scalarser(gen_1, 0, rk): zeroser(0, 0);
421 63 : setvalser(t, 1);
422 63 : Lrk = ellL1_der(e, vec, &C, t, rk, prec);
423 63 : if (DEBUGLEVEL) timer_printf(&ti, "L^(%ld)=%Ps", rk, Lrk);
424 63 : if (abscmprr(Lrk, eps) > 0)
425 35 : return gc_GEN(av, mkvec2(stoi(rk), Lrk));
426 28 : set_avma(av2);
427 : }
428 : }
429 :
430 : /* Heegner point computation
431 :
432 : This section is a C version by Bill Allombert of a GP script by
433 : Christophe Delaunay which was based on a GP script by John Cremona.
434 : Reference: Henri Cohen's book GTM 239.
435 : */
436 :
437 : static void
438 0 : heegner_L1_bg(void*E, GEN n, GEN a)
439 : {
440 0 : struct baby_giant *bb = (struct baby_giant*) E;
441 0 : long j, l = lg(bb->giant);
442 0 : for (j = 1; j < l; j++)
443 0 : if (cmpii(n, gel(bb->bnd,j)) <= 0)
444 : {
445 0 : ulong r, q = uabsdiviu_rem(n, bb->rbnd[j], &r);
446 0 : GEN giant = gel(bb->giant, j), baby = gel(bb->baby, j);
447 0 : affgc(gadd(gel(giant, q+1), gdiv(gmul(gel(baby, r+1), a), n)),
448 0 : gel(giant, q+1));
449 : }
450 0 : }
451 :
452 : static void
453 2869580 : heegner_L1(void*E, GEN n, GEN a)
454 : {
455 2869580 : struct baby_giant *bb = (struct baby_giant*) E;
456 2869580 : long j, l = lg(bb->giant);
457 16808722 : for (j = 1; j < l; j++)
458 13939142 : if (cmpii(n, gel(bb->bnd,j)) <= 0)
459 : {
460 11892741 : ulong r, q = uabsdiviu_rem(n, bb->rbnd[j], &r);
461 11892741 : GEN giant = gel(bb->giant, j), baby = gel(bb->baby, j);
462 11892741 : GEN ex = mulreal(gel(baby, r+1), gel(giant, q+1));
463 11892741 : affrr(addrr(gel(bb->sum, j), divri(mulri(ex, a), n)),
464 11892741 : gel(bb->sum, j));
465 : }
466 2869580 : }
467 : /* export ? */
468 : static GEN
469 0 : ctoc(GEN x, long prec)
470 0 : { GEN y = cgetc(prec); affgc(x, y); return y; }
471 :
472 : /* [powers(x[i], n[i]), i=1..#x] */
473 : static GEN
474 42 : RgV_zv_powers(GEN x, GEN n)
475 189 : { pari_APPLY_same(gpowers(gel(x,i), n[i])); }
476 :
477 : /* Return C, C[i][j] = Q[j]^i, i = 1..nb */
478 : static void
479 0 : baby_init2(struct baby_giant *bb, GEN Q, GEN bnd, GEN rbnd, long prec)
480 : {
481 0 : long i, j, l = lg(Q);
482 0 : GEN C = RgV_zv_powers(Q, rbnd), R = cgetg(l,t_VEC);
483 0 : for (i = 1; i < l; ++i)
484 : {
485 0 : gel(R, i) = cgetg(rbnd[i]+1, t_VEC);
486 0 : gmael(R, i, 1) = ctoc(gmael(C, i, 2), prec);
487 0 : for (j = 2; j <= rbnd[i]; j++) gmael(R, i, j) = ctoc(gen_0, prec);
488 : }
489 0 : bb->baby = C; bb->giant = R;
490 0 : bb->bnd = bnd; bb->rbnd = rbnd;
491 0 : }
492 :
493 : /* Return C, C[i][j] = Q[j]^i, i = 1..nb */
494 : static void
495 42 : baby_init3(struct baby_giant *bb, GEN Q, GEN bnd, GEN rbnd, long prec)
496 : {
497 42 : long i, l = lg(Q);
498 42 : GEN S, C = RgV_zv_powers(Q, rbnd), R = cgetg(l,t_VEC);
499 189 : for (i = 1; i < l; ++i)
500 147 : gel(R, i) = gpowers(gmael(C, i, 1+rbnd[i]), rbnd[i]);
501 42 : S = cgetg(l,t_VEC);
502 189 : for (i = 1; i < l; ++i) gel(S, i) = rtor(real_i(gmael(C, i, 2)), prec);
503 42 : bb->baby = C; bb->giant = R; bb->sum = S;
504 42 : bb->bnd = bnd; bb->rbnd = rbnd;
505 42 : }
506 :
507 : /* ymin a t_REAL */
508 : static GEN
509 42 : heegner_psi(GEN E, GEN N, GEN points, long bitprec)
510 : {
511 42 : pari_sp av = avma, av2;
512 : struct baby_giant bb;
513 : struct bg_data bg;
514 42 : long l, k, L = lg(points)-1, prec = nbits2prec(bitprec)+EXTRAPREC64;
515 42 : GEN Q, pi2 = Pi2n(1, prec), bnd, rbnd, bndmax;
516 42 : GEN B = divrr(mulur(bitprec,mplog2(DEFAULTPREC)), pi2);
517 :
518 42 : rbnd = cgetg(L+1, t_VECSMALL); av2 = avma;
519 42 : bnd = cgetg(L+1, t_VEC);
520 42 : Q = cgetg(L+1, t_VEC);
521 189 : for (l = 1; l <= L; ++l)
522 : {
523 147 : gel(bnd,l) = ceil_safe(divrr(B,imag_i(gel(points, l))));
524 147 : rbnd[l] = itou(sqrtint(gel(bnd,l)))+1;
525 147 : gel(Q, l) = expIxy(pi2, gel(points, l), prec);
526 : }
527 42 : (void)gc_all(av2, 2, &bnd, &Q);
528 42 : bndmax = gel(bnd,vecindexmax(bnd));
529 42 : gen_BG_init(&bg, E, N, bndmax);
530 42 : if (bitprec >= 1900)
531 : {
532 0 : GEN S = cgetg(L+1, t_VEC);
533 0 : baby_init2(&bb, Q, bnd, rbnd, prec);
534 0 : gen_BG_rec((void*)&bb, heegner_L1_bg, &bg);
535 0 : for (k = 1; k <= L; ++k)
536 : {
537 0 : pari_sp av2 = avma;
538 0 : long j, g = rbnd[k];
539 0 : GEN giant = gmael(bb.baby, k, g+1), Sl = real_0(prec);
540 0 : for (j = g; j >=1; j--) Sl = gadd(gmul(Sl, giant), gmael(bb.giant,k,j));
541 0 : gel(S, k) = gc_upto(av2, real_i(Sl));
542 : }
543 0 : return gc_upto(av, S);
544 : }
545 : else
546 : {
547 42 : baby_init3(&bb, Q, bnd, rbnd, prec);
548 42 : gen_BG_rec((void*)&bb, heegner_L1, &bg);
549 42 : return gc_GEN(av, bb.sum);
550 : }
551 : }
552 :
553 : /*Returns lambda_bad list for one prime p, nv = localred(E, p) */
554 : static GEN
555 91 : lambda1(GEN E, GEN nv, GEN p, long prec)
556 : {
557 : pari_sp av;
558 : GEN res, lp;
559 91 : long kod = itos(gel(nv, 2));
560 91 : if (kod==2 || kod ==-2) return cgetg(1,t_VEC);
561 91 : av = avma; lp = glog(p, prec);
562 91 : if (kod > 4)
563 : {
564 14 : long n = Z_pval(ell_get_disc(E), p);
565 14 : long j, m = kod - 4, nl = 1 + (m >> 1L);
566 14 : res = cgetg(nl, t_VEC);
567 35 : for (j = 1; j < nl; j++)
568 21 : gel(res, j) = gmul(lp, gsubgs(gdivgu(sqru(j), n), j)); /* j^2/n - j */
569 : }
570 77 : else if (kod < -4)
571 14 : res = mkvec2(negr(lp), shiftr(mulrs(lp, kod), -2));
572 : else
573 : {
574 63 : const long lam[] = {8,9,0,6,0,0,0,3,4};
575 63 : long m = -lam[kod+4];
576 63 : res = mkvec(divru(mulrs(lp, m), 6));
577 : }
578 91 : return gc_GEN(av, res);
579 : }
580 :
581 : static GEN
582 42 : lambdalist(GEN E, long prec)
583 : {
584 42 : pari_sp ltop = avma;
585 42 : GEN glob = ellglobalred(E), plist = gmael(glob,4,1), L = gel(glob,5);
586 42 : GEN res, v, D = ell_get_disc(E);
587 42 : long i, j, k, l, m, n, np = lg(plist), lr = 1;
588 42 : v = cgetg(np, t_VEC);
589 147 : for (j = 1, i = 1 ; j < np; ++j)
590 : {
591 105 : GEN p = gel(plist, j);
592 105 : if (dvdii(D, sqri(p)))
593 : {
594 91 : GEN la = lambda1(E, gel(L,j), p, prec);
595 91 : gel(v, i++) = la;
596 91 : lr *= lg(la);
597 : }
598 : }
599 42 : np = i;
600 42 : res = cgetg(lr+1, t_VEC);
601 42 : gel(res, 1) = gen_0; n = 1; m = 1;
602 133 : for (j = 1; j < np; ++j)
603 : {
604 91 : GEN w = gel(v, j);
605 91 : long lw = lg(w);
606 308 : for (k = 1; k <= n; k++)
607 : {
608 217 : GEN t = gel(res, k);
609 455 : for (l = 1, m = n; l < lw; l++, m+=n)
610 238 : gel(res, k + m) = mpadd(t, gel(w, l));
611 : }
612 91 : n = m;
613 : }
614 42 : return gc_GEN(ltop, res);
615 : }
616 :
617 : /* P a t_INT or t_FRAC, return its logarithmic height */
618 : static GEN
619 98 : heightQ(GEN P, long prec)
620 : {
621 : long s;
622 98 : if (typ(P) == t_FRAC)
623 : {
624 56 : GEN a = gel(P,1), b = gel(P,2);
625 56 : P = abscmpii(a,b) > 0 ? a: b;
626 : }
627 98 : s = signe(P);
628 98 : if (!s) return real_0(prec);
629 84 : if (s < 0) P = negi(P);
630 84 : return glog(P, prec);
631 : }
632 :
633 : /* t a t_INT or t_FRAC, returns max(1, log |t|), returns a t_REAL */
634 : static GEN
635 147 : logplusQ(GEN t, long prec)
636 : {
637 147 : if (typ(t) == t_INT)
638 : {
639 42 : if (!signe(t)) return real_1(prec);
640 28 : t = absi_shallow(t);
641 : }
642 : else
643 : {
644 105 : GEN a = gel(t,1), b = gel(t,2);
645 105 : if (abscmpii(a, b) < 0) return real_1(prec);
646 56 : if (signe(a) < 0) t = gneg(t);
647 : }
648 84 : return glog(t, prec);
649 : }
650 :
651 : /* See GTM239, p532, Th 8.1.18
652 : * Return M such that h_naive <= M */
653 : GEN
654 98 : hnaive_max(GEN ell, GEN ht)
655 : {
656 98 : const long prec = LOWDEFAULTPREC; /* minimal accuracy */
657 98 : GEN b2 = ell_get_b2(ell), j = ell_get_j(ell);
658 98 : GEN logd = glog(absi_shallow(ell_get_disc(ell)), prec);
659 98 : GEN logj = logplusQ(j, prec);
660 98 : GEN hj = heightQ(j, prec);
661 49 : GEN logb2p = signe(b2)? addrr(logplusQ(gdivgu(b2, 12),prec), mplog2(prec))
662 98 : : real_1(prec);
663 98 : GEN mu = addrr(divru(addrr(logd, logj),6), logb2p);
664 98 : return addrs(addrr(addrr(ht, divru(hj,12)), mu), 2);
665 : }
666 :
667 : static GEN
668 147 : qfb_root(GEN Q, GEN vDi)
669 : {
670 147 : GEN a2 = shifti(gel(Q, 1),1), b = gel(Q, 2);
671 147 : return mkcomplex(gdiv(negi(b),a2),divri(vDi,a2));
672 : }
673 :
674 : static GEN
675 24668 : qimag2(GEN Q)
676 : {
677 24668 : pari_sp av = avma;
678 24668 : GEN z = gdiv(negi(qfb_disc(Q)), shifti(sqri(gel(Q, 1)),2));
679 24668 : return gc_upto(av, z);
680 : }
681 :
682 : /***************************************************/
683 : /*Routines for increasing the imaginary parts using*/
684 : /*Atkin-Lehner operators */
685 : /***************************************************/
686 :
687 : static GEN
688 24668 : qfb_mult(GEN Q, GEN a, GEN b, GEN c, GEN d)
689 : {
690 24668 : GEN A = gel(Q, 1) , B = gel(Q, 2), C = gel(Q, 3), D = qfb_disc(Q);
691 24668 : GEN a2 = sqri(a), b2 = sqri(b), c2 = sqri(c), d2 = sqri(d);
692 24668 : GEN ad = mulii(d, a), bc = mulii(b, c), e = subii(ad, bc);
693 24668 : GEN W1 = addii(addii(mulii(a2, A), mulii(mulii(c, a), B)), mulii(c2, C));
694 24668 : GEN W3 = addii(addii(mulii(b2, A), mulii(mulii(d, b), B)), mulii(d2, C));
695 24668 : GEN W2 = addii(addii(mulii(mulii(shifti(b,1), a), A),
696 : mulii(addii(ad, bc), B)),
697 : mulii(mulii(shifti(d,1), c), C));
698 24668 : if (!equali1(e)) {
699 22190 : W1 = diviiexact(W1,e);
700 22190 : W2 = diviiexact(W2,e);
701 22190 : W3 = diviiexact(W3,e);
702 : }
703 24668 : return mkqfb(W1, W2, W3, D);
704 : }
705 :
706 : #ifdef DEBUG
707 : static void
708 : best_point_old(GEN Q, GEN NQ, GEN f, GEN *u, GEN *v)
709 : {
710 : long n, k;
711 : GEN U, c, d, A = gel(f,1), B = gel(f,2), C = gel(f,3), D = qfb_disc(f);
712 : GEN q = mkqfb(mulii(NQ, C), negi(B), diviiexact(A, NQ), D);
713 : redimagsl2(q, &U);
714 : *u = c = gcoeff(U, 1, 1);
715 : *v = d = gcoeff(U, 2, 1);
716 : if (equali1(gcdii(mulii(*u, NQ), mulii(*v, Q)))) return;
717 : for (n = 1;; n++)
718 : {
719 : for (k = -n; k <= n; k++)
720 : {
721 : *u = addis(c, k); *v = addiu(d, n);
722 : if (equali1(gcdii(mulii(*u, NQ), mulii(*v, Q)))) return;
723 : *v = subiu(d, n);
724 : if (equali1(gcdii(mulii(*u, NQ), mulii(*v, Q)))) return;
725 : *u = addiu(c, n); *v = addis(d, k);
726 : if (equali1(gcdii(mulii(*u, NQ), mulii(*v, Q)))) return;
727 : *u = subiu(c, n);
728 : if (equali1(gcdii(mulii(*u, NQ), mulii(*v, Q)))) return;
729 : }
730 : }
731 : }
732 : /* q(x,y) = ax^2 + bxy + cy^2 */
733 : static GEN
734 : qfb_eval(GEN q, GEN x, GEN y)
735 : {
736 : GEN a = gel(q,1), b = gel(q,2), c = gel(q,3);
737 : GEN x2 = sqri(x), y2 = sqri(y), xy = mulii(x,y);
738 : return addii(addii(mulii(a, x2), mulii(b,xy)), mulii(c, y2));
739 : }
740 : #endif
741 :
742 : static long
743 6580 : nexti(long i) { return i>0 ? -i : 1-i; }
744 :
745 : /* q0 + i q1 + i^2 q2 */
746 : static GEN
747 12313 : qfmin_eval(GEN q0, GEN q1, GEN q2, long i)
748 12313 : { return addii(mulis(addii(mulis(q2, i), q1), i), q0); }
749 :
750 : /* assume a > 0, return gcd(a,b,c) */
751 : static ulong
752 16436 : gcduii(ulong a, GEN b, GEN c)
753 : {
754 16436 : a = ugcdiu(b, a);
755 16436 : return a == 1? 1: ugcdiu(c, a);
756 : }
757 :
758 : static void
759 24668 : best_point(GEN Q, GEN NQ, GEN f, GEN *pu, GEN *pv)
760 : {
761 24668 : GEN a = mulii(NQ, gel(f,3)), b = negi(gel(f,2)), c = diviiexact(gel(f,1), NQ);
762 24668 : GEN D = qfb_disc(f);
763 24668 : GEN U, qr = redimagsl2(mkqfb(a, b, c, D), &U);
764 24668 : GEN A = gel(qr,1), B = gel(qr,2), A2 = shifti(A,1), AA4 = sqri(A2);
765 : GEN V, best;
766 : long y;
767 :
768 24668 : D = absi_shallow(D);
769 : /* 4A qr(x,y) = (2A x + By)^2 + D y^2
770 : * Write x = x0(y) + i, where x0 is an integer minimum
771 : * (the smallest in case of tie) of x-> qr(x,y), for given y.
772 : * 4A qr(x,y) = ((2A x0 + By)^2 + Dy^2) + 4A i (2A x0 + By) + 4A^2 i^2
773 : * = q0(y) + q1(y) i + q2 i^2
774 : * Loop through (x,y), y>0 by (roughly) increasing values of qr(x,y) */
775 :
776 : /* We must test whether [X,Y]~ := U * [x,y]~ satisfy (X NQ, Y Q) = 1
777 : * This is equivalent to (X,Y) = 1 (note that (X,Y) = (x,y)), and
778 : * (X, Q) = (Y, NQ) = 1.
779 : * We have U * [x0+i, y]~ = U * [x0,y]~ + i U[,1] =: V0 + i U[,1] */
780 :
781 : /* try [1,0]~ = first minimum */
782 24668 : V = gel(U,1); /* U *[1,0]~ */
783 24668 : *pu = gel(V,1);
784 24668 : *pv = gel(V,2);
785 30947 : if (is_pm1(gcdii(*pu, Q)) && is_pm1(gcdii(*pv, NQ))) return;
786 :
787 : /* try [0,1]~ = second minimum */
788 11935 : V = gel(U,2); /* U *[0,1]~ */
789 11935 : *pu = gel(V,1);
790 11935 : *pv = gel(V,2);
791 11935 : if (is_pm1(gcdii(*pu, Q)) && is_pm1(gcdii(*pv, NQ))) return;
792 :
793 : /* (X,Y) = (1, \pm1) always works. Try to do better now */
794 5656 : best = subii(addii(a, c), absi_shallow(b));
795 5656 : *pu = gen_1;
796 5656 : *pv = signe(b) < 0? gen_1: gen_m1;
797 :
798 5656 : for (y = 1;; y++)
799 9128 : {
800 : GEN Dy2, r, By, x0, q0, q1, V0;
801 : long i;
802 14784 : if (y > 1)
803 : {
804 10703 : if (gcduii(y, gcoeff(U,1,1), Q) != 1) continue;
805 7308 : if (gcduii(y, gcoeff(U,2,1), NQ) != 1) continue;
806 : }
807 11396 : Dy2 = mulii(D, sqru(y));
808 11396 : if (cmpii(Dy2, best) >= 0) break; /* we won't improve. STOP */
809 5740 : By = muliu(B,y), x0 = truedvmdii(negi(By), A2, &r);
810 5740 : if (cmpii(r, A) >= 0) { x0 = subiu(x0,1); r = subii(r, A2); }
811 : /* (2A x + By)^2 + Dy^2, minimal at x = x0. Assume A2 > 0 */
812 : /* r = 2A x0 + By */
813 5740 : q0 = addii(sqri(r), Dy2); /* minimal value for this y, at x0 */
814 5740 : if (cmpii(q0, best) >= 0) continue; /* we won't improve for this y */
815 5733 : q1 = shifti(mulii(A2, r), 1);
816 :
817 5733 : V0 = ZM_ZC_mul(U, mkcol2(x0, utoipos(y)));
818 12313 : for (i = 0;; i = nexti(i))
819 6580 : {
820 12313 : pari_sp av2 = avma;
821 12313 : GEN x, N = qfmin_eval(q0, q1, AA4, i);
822 12313 : if (cmpii(N , best) >= 0) break;
823 12271 : x = addis(x0, i);
824 12271 : if (ugcdiu(x, y) == 1)
825 : {
826 : GEN u, v;
827 12229 : V = ZC_add(V0, ZC_z_mul(gel(U,1), i)); /* [X, Y] */
828 12229 : u = gel(V,1);
829 12229 : v = gel(V,2);
830 12229 : if (is_pm1(gcdii(u, Q)) && is_pm1(gcdii(v, NQ)))
831 : {
832 5691 : *pu = u;
833 5691 : *pv = v;
834 5691 : best = N; break;
835 : }
836 : }
837 6580 : set_avma(av2);
838 : }
839 : }
840 : #ifdef DEBUG
841 : {
842 : GEN oldu, oldv, F = mkqfb(a, b, c, qfb_disc(f));
843 : best_point_old(Q, NQ, f, &oldu, &oldv);
844 : if (!equalii(oldu, *pu) || !equalii(oldv, *pv))
845 : {
846 : if (!equali1(gcdii(mulii(*pu, NQ), mulii(*pv, Q))))
847 : pari_err_BUG("best_point (gcd)");
848 : if (cmpii(qfb_eval(F, *pu,*pv), qfb_eval(F, oldu, oldv)) > 0)
849 : {
850 : pari_warn(warner, "%Ps,%Ps,%Ps, %Ps > %Ps",
851 : Q,NQ,f, mkvec2(*pu,*pv), mkvec2(oldu,oldv));
852 : pari_err_BUG("best_point (too large)");
853 : }
854 : }
855 : }
856 : #endif
857 : }
858 :
859 : static GEN
860 24668 : best_lift(GEN Q, GEN NQ, GEN f)
861 : {
862 : GEN a, b, c, d, dQ, cNQ;
863 24668 : best_point(Q, NQ, f, &c, &d);
864 24668 : dQ = mulii(d, Q); cNQ = mulii(NQ, c);
865 24668 : (void)bezout(dQ, cNQ, &a, &b);
866 24668 : return qfb_mult(f, dQ, b, mulii(negi(Q),cNQ), mulii(a,Q));
867 : }
868 :
869 : static GEN
870 2478 : lift_points(GEN listQ, GEN f, GEN *pt, GEN *pQ)
871 : {
872 2478 : pari_sp av = avma;
873 2478 : GEN yf = gen_0, tf = NULL, Qf = NULL;
874 2478 : long k, l = lg(listQ);
875 27146 : for (k = 1; k < l; ++k)
876 : {
877 24668 : GEN c = gel(listQ, k), Q = gel(c,1), NQ = gel(c,2);
878 24668 : GEN t = best_lift(Q, NQ, f), y = qimag2(t);
879 24668 : if (gcmp(y, yf) > 0) { yf = y; Qf = Q; tf = t; }
880 : }
881 2478 : *pt = tf; *pQ = Qf; return gc_all(av, 3, &yf, pt, pQ);
882 : }
883 :
884 : /***************************/
885 : /* Twists */
886 : /***************************/
887 :
888 : static GEN
889 56 : ltwist1(GEN E, GEN d, long bitprec)
890 : {
891 56 : pari_sp av = avma;
892 56 : GEN Ed = elltwist(E, d), z = ellL1(Ed, 0, bitprec);
893 56 : obj_free(Ed); return gc_uptoleaf(av, z);
894 : }
895 :
896 : /* Return O_re*c(E)/(4*O_vol*|E_t|^2) */
897 :
898 : static GEN
899 42 : heegner_indexmult(GEN om, long t, GEN tam, long prec)
900 : {
901 42 : pari_sp av = avma;
902 42 : GEN Ovr = gabs(imag_i(gel(om, 2)), prec); /* O_vol/O_re, t_REAL */
903 42 : return gc_upto(av, divru(divir(tam, Ovr), 4*t*t));
904 : }
905 :
906 : /* omega(gcd(D, N)), given faN = factor(N) */
907 : static long
908 56 : omega_N_D(GEN faN, ulong D)
909 : {
910 56 : GEN P = gel(faN, 1);
911 56 : long i, l = lg(P), w = 0;
912 196 : for (i = 1; i < l; i++)
913 140 : if (dvdui(D, gel(P,i))) w++;
914 56 : return w;
915 : }
916 :
917 : static GEN
918 56 : heegner_indexmultD(GEN faN, GEN a, long D, GEN sqrtD)
919 : {
920 56 : pari_sp av = avma;
921 : GEN c;
922 : long w;
923 56 : switch(D)
924 : {
925 0 : case -3: w = 9; break;
926 0 : case -4: w = 4; break;
927 56 : default: w = 1;
928 : }
929 56 : c = shifti(stoi(w), omega_N_D(faN,-D)); /* (w(D)/2)^2 * 2^omega(gcd(D,N)) */
930 56 : return gc_upto(av, mulri(mulrr(a, sqrtD), c));
931 : }
932 :
933 : static GEN
934 399 : nf_to_basis(GEN nf, GEN x)
935 : {
936 399 : x = nf_to_scalar_or_basis(nf, x);
937 399 : if (typ(x)!=t_COL)
938 301 : x = scalarcol(x, nf_get_degree(nf));
939 399 : return x;
940 : }
941 :
942 : static GEN
943 196 : etnf_to_basis(GEN et, GEN x)
944 : {
945 196 : long i, l = lg(et);
946 196 : GEN V = cgetg(l, t_VEC);
947 595 : for (i = 1; i < l; i++)
948 399 : gel(V,i) = nf_to_basis(gel(et,i), x);
949 196 : return shallowconcat1(V);
950 : }
951 :
952 : static GEN
953 140 : etnf_get_M(GEN et)
954 : {
955 140 : long i, l = lg(et);
956 140 : GEN V = cgetg(l, t_VEC);
957 448 : for (i = 1; i < l; i++)
958 308 : gel(V,i)=nf_get_M(gel(et,i));
959 140 : return shallowmatconcat(diagonal(V));
960 : }
961 :
962 : static long
963 49 : etnf_get_varn(GEN et)
964 : {
965 49 : return nf_get_varn(gel(et,1));
966 : }
967 :
968 : static GEN
969 98 : heegner_descent_try_point(GEN nfA, GEN z, GEN den, long prec)
970 : {
971 98 : pari_sp av = avma;
972 98 : GEN etal = gel(nfA,1), A = gel(nfA,2), cb = gel(nfA,3);
973 98 : GEN al = gel(nfA,4), th = gel(nfA, 5);
974 98 : GEN et = gel(etal,1), zk = gel(etal, 2), T = gel(etal,3);
975 98 : GEN M = etnf_get_M(et);
976 98 : long i, j, n = lg(th)-1, l = lg(al);
977 98 : GEN u2 = gsqr(gel(cb,1)), r = gel(cb,2);
978 98 : GEN zz = gdiv(gsub(z,r), u2);
979 98 : GEN be = cgetg(n+1, t_COL);
980 154 : for (j = 1; j < l; j++)
981 : {
982 98 : GEN aj = gel(al, j), Aj = gel(A,j);
983 371 : for (i = 1; i <= n; i++)
984 273 : gel(be,i) = gsqrt(gmul(gsub(zz, gel(th,i)), gel(aj,i)), prec);
985 336 : for (i = 0; i <= (1<<(n-1))-1; i++)
986 : {
987 : long eps;
988 280 : GEN s = gmul(den, RgM_solve_realimag(M, be));
989 280 : GEN S = grndtoi(s, &eps), V, S2;
990 280 : gel(be,1+odd(i)) = gneg(gel(be,1+odd(i)));
991 280 : if (eps > -7) continue;
992 42 : S2 = QXQ_sqr(RgV_RgC_mul(zk, S), T);
993 42 : V = gdiv(QXQ_mul(S2, Aj, T), sqri(den));
994 42 : if (typ(V) != t_POL || degpol(V) != 1) continue;
995 42 : if (gequalm1(gel(V,3)))
996 42 : return gc_upto(av,gadd(gmul(gel(V,2),u2),r));
997 : }
998 : }
999 56 : return gc_NULL(av);
1000 : }
1001 :
1002 : static GEN
1003 1785 : heegner_try_point(GEN E, GEN nfA, GEN lambdas, GEN ht, GEN z, long prec)
1004 : {
1005 1785 : long l = lg(lambdas);
1006 : long i, eps;
1007 1785 : GEN P = real_i(pointell(E, z, prec)), x = gel(P,1);
1008 1785 : GEN rh = subrr(ht, shiftr(ellheightoo(E, P, prec),1));
1009 26572 : for (i = 1; i < l; ++i)
1010 : {
1011 24829 : GEN logd = shiftr(gsub(rh, gel(lambdas, i)), -1);
1012 24829 : GEN d, approxd = gexp(logd, prec);
1013 24829 : d = grndtoi(approxd, &eps);
1014 24829 : if (signe(d) > 0 && eps<-10)
1015 : {
1016 : GEN X, ylist;
1017 98 : if (DEBUGLEVEL > 2)
1018 0 : err_printf("\nTrying lambda number %ld, logd=%Ps, approxd=%Ps\n", i, logd, approxd);
1019 98 : X = heegner_descent_try_point(nfA, x, d, prec);
1020 98 : if (X)
1021 : {
1022 42 : ylist = ellordinate(E, X, prec);
1023 42 : if (lg(ylist) > 1)
1024 : {
1025 42 : GEN P = mkvec2(X, gel(ylist, 1));
1026 42 : GEN hp = ellheight(E,P,prec);
1027 42 : if (signe(hp) && cmprr(hp, shiftr(ht,1)) < 0 && cmprr(hp, shiftr(ht,-1)) > 0)
1028 42 : return P;
1029 0 : if (DEBUGLEVEL)
1030 0 : err_printf("found non-Heegner point %Ps\n", P);
1031 : }
1032 : }
1033 : }
1034 : }
1035 1743 : return NULL;
1036 : }
1037 :
1038 : static GEN
1039 42 : heegner_find_point(GEN e, GEN nfA, GEN om, GEN ht, GEN z1, long k, long prec)
1040 : {
1041 42 : GEN lambdas = lambdalist(e, prec);
1042 42 : pari_sp av = avma;
1043 : long m;
1044 42 : GEN Ore = gel(om, 1), Oim = gel(om, 2);
1045 42 : if (DEBUGLEVEL)
1046 0 : err_printf("%ld*%ld multipliers to test: ",k,lg(lambdas)-1);
1047 966 : for (m = 0; m < k; m++)
1048 : {
1049 966 : GEN P, z2 = divru(addrr(z1, mulsr(m, Ore)), k);
1050 966 : if (DEBUGLEVEL > 2)
1051 0 : err_printf("%ld ",m);
1052 966 : P = heegner_try_point(e, nfA, lambdas, ht, z2, prec);
1053 966 : if (P) return P;
1054 931 : if (signe(ell_get_disc(e)) > 0)
1055 : {
1056 819 : z2 = gadd(z2, gmul2n(Oim, -1));
1057 819 : P = heegner_try_point(e, nfA, lambdas, ht, z2, prec);
1058 819 : if (P) return P;
1059 : }
1060 924 : set_avma(av);
1061 : }
1062 0 : pari_err_BUG("ellheegner, point not found");
1063 : return NULL; /* LCOV_EXCL_LINE */
1064 : }
1065 :
1066 : /* N > 1, fa = factor(N), return factor(4*N) */
1067 : static GEN
1068 42 : fa_shift2(GEN fa)
1069 : {
1070 42 : GEN P = gel(fa,1), E = gel(fa,2);
1071 42 : if (absequaliu(gcoeff(fa,1,1), 2))
1072 : {
1073 21 : E = shallowcopy(E);
1074 21 : gel(E,1) = addiu(gel(E,1), 2);
1075 : }
1076 : else
1077 : {
1078 21 : P = shallowconcat(gen_2, P);
1079 21 : E = shallowconcat(gen_2, E);
1080 : }
1081 42 : return mkmat2(P, E);
1082 : }
1083 :
1084 : /* P = prime divisors of N(E). Return the product of primes p in P, a_p != -1
1085 : * HACK: restrict to small primes since large ones won't divide our C-long
1086 : * discriminants */
1087 : static GEN
1088 42 : get_bad(GEN E, GEN P)
1089 : {
1090 42 : long k, l = lg(P), ibad = 1;
1091 42 : GEN B = cgetg(l, t_VECSMALL);
1092 147 : for (k = 1; k < l; k++)
1093 : {
1094 105 : GEN p = gel(P,k);
1095 105 : long pp = itos_or_0(p);
1096 105 : if (!pp) break;
1097 105 : if (! equalim1(ellap(E,p))) B[ibad++] = pp;
1098 : }
1099 42 : setlg(B, ibad); return ibad == 1? NULL: zv_prod_Z(B);
1100 : }
1101 :
1102 : /* list of pairs [Q,N/Q], where Q | N and gcd(Q,N/Q) = 1 */
1103 : static GEN
1104 42 : find_div(GEN N, GEN faN)
1105 : {
1106 42 : GEN listQ = divisors(faN);
1107 42 : long j, k, l = lg(listQ);
1108 :
1109 42 : gel(listQ, 1) = mkvec2(gen_1, N);
1110 1624 : for (j = k = 2; k < l; ++k)
1111 : {
1112 1582 : GEN Q = gel(listQ, k), NQ = diviiexact(N, Q);
1113 1582 : if (is_pm1(gcdii(Q,NQ))) gel(listQ, j++) = mkvec2(Q,NQ);
1114 : }
1115 42 : setlg(listQ, j); return listQ;
1116 : }
1117 :
1118 : static long
1119 8652 : testDisc(GEN bad, long d) { return !bad || ugcdiu(bad, -d) == 1; }
1120 : /* bad = product of bad primes. Return the NDISC largest fundamental
1121 : * discriminants D < d such that (D,bad) = 1 and d is a square mod 4N */
1122 : static GEN
1123 42 : listDisc(GEN fa4N, GEN bad, long d, long ndisc)
1124 : {
1125 42 : GEN v = cgetg(ndisc+1, t_VECSMALL);
1126 42 : pari_sp av = avma;
1127 42 : long j = 1;
1128 : for(;;)
1129 : {
1130 8652 : d -= odd(d)? 1: 3;
1131 8652 : if (testDisc(bad,d) && unegisfundamental(-d) && Zn_issquare(stoi(d), fa4N))
1132 : {
1133 420 : v[j++] = d;
1134 420 : if (j > ndisc) break;
1135 : }
1136 8610 : set_avma(av);
1137 : }
1138 42 : set_avma(av); return v;
1139 : }
1140 : /* L = vector of [q1,q2] or [q1,q2,q2']
1141 : * cd = (b^2 - D)/(4N) */
1142 : static void
1143 166880 : listfill(GEN N, GEN b, GEN c, GEN d, GEN D, GEN L, long *s)
1144 : {
1145 166880 : long k, l = lg(L);
1146 166880 : GEN add, frm2, a = mulii(d, N), V = mkqfb(a,b,c,D), frm = qfbred_i(V);
1147 600089 : for (k = 1; k < l; ++k)
1148 : { /* Lk = [v,frm] or [v,frm,frm2] */
1149 597611 : GEN Lk = gel(L,k);
1150 : long i;
1151 1509515 : for (i = 2; i < lg(Lk); i++) /* 1 or 2 elements */
1152 1076306 : if (gequal(frm, gel(Lk,i)))
1153 : {
1154 164402 : GEN v = gel(Lk, 1);
1155 164402 : if (cmpii(a, gel(v,1)) < 0) gel(Lk,1) = V;
1156 164402 : return;
1157 : }
1158 : }
1159 2478 : frm2 = qfbred_i(mkqfb(d, negi(b), mulii(c,N), D));
1160 2478 : add = gequal(frm, frm2)? mkvec2(V,frm): mkvec3(V,frm,frm2);
1161 2478 : vectrunc_append(L, add);
1162 2478 : *s += lg(add) - 2;
1163 : }
1164 : /* faN4 = factor(4*N) */
1165 : static GEN
1166 420 : listheegner(GEN N, GEN faN4, GEN listQ, GEN D)
1167 : {
1168 420 : pari_sp av = avma;
1169 420 : const long kmin = 30;
1170 420 : long h = itos(quadclassno(D));
1171 420 : GEN ymin, b = Zn_sqrt(D, faN4), L = vectrunc_init(h+1);
1172 420 : long l, k, s = 0;
1173 13020 : for (k = 0; k < kmin || s < h; k++)
1174 : {
1175 12600 : GEN bk = addii(b, mulsi(2*k, N));
1176 12600 : GEN C = diviiexact(shifti(subii(sqri(bk), D), -2), N);
1177 12600 : GEN div = divisors(C);
1178 12600 : long i, l = lg(div);
1179 179480 : for (i = 1; i < l; i++)
1180 : {
1181 166880 : GEN d = gel(div, i), c = gel(div, l-i); /* cd = C */
1182 166880 : listfill(N, bk, c, d, D, L, &s);
1183 : }
1184 : }
1185 420 : l = lg(L); ymin = NULL;
1186 2898 : for (k = 1; k < l; k++)
1187 : {
1188 2478 : GEN t, Q, Lk = gel(L,k), f = gel(Lk,1);
1189 2478 : GEN y = lift_points(listQ, f, &t, &Q);
1190 2478 : gel(L, k) = mkvec3(t, stoi(lg(Lk) - 2), Q);
1191 2478 : if (!ymin || gcmp(y, ymin) < 0) ymin = y;
1192 : }
1193 420 : if (DEBUGLEVEL > 1)
1194 0 : err_printf("Disc %Ps : N*ymin = %Pg\n", D,
1195 : gmul(gsqrt(ymin, DEFAULTPREC),N));
1196 420 : return gc_GEN(av, mkvec3(ymin, L, D));
1197 : }
1198 :
1199 : /* Q | N, P = prime divisors of N, R[i] = local epsilon-factor at P[i].
1200 : * Return \prod_{p | Q} R[i] */
1201 : static long
1202 147 : rootno(GEN Q, GEN P, GEN R)
1203 : {
1204 147 : long s = 1, i, l = lg(P);
1205 581 : for (i = 1; i < l; i++)
1206 434 : if (dvdii(Q, gel(P,i))) s *= R[i];
1207 147 : return s;
1208 : }
1209 :
1210 : static void
1211 42 : heegner_find_disc(GEN *points, GEN *coefs, long *pind, GEN E,
1212 : GEN indmult, long ndisc, long prec)
1213 : {
1214 42 : long d = 0;
1215 : GEN faN4, bad, N, faN, listQ, listR;
1216 :
1217 42 : ellQ_get_Nfa(E, &N, &faN);
1218 42 : faN4 = fa_shift2(faN);
1219 42 : listQ = find_div(N, faN);
1220 42 : bad = get_bad(E, gel(faN, 1));
1221 42 : listR = gel(obj_check(E, Q_ROOTNO), 2);
1222 : for(;;)
1223 0 : {
1224 42 : pari_sp av = avma;
1225 42 : GEN list, listD = listDisc(faN4, bad, d, ndisc);
1226 42 : long k, l = lg(listD);
1227 42 : list = cgetg(l, t_VEC);
1228 462 : for (k = 1; k < l; ++k)
1229 420 : gel(list, k) = listheegner(N, faN4, listQ, stoi(listD[k]));
1230 42 : list = vecsort0(list, gen_1, 0);
1231 56 : for (k = l-1; k > 0; --k)
1232 : {
1233 56 : long bprec = 8;
1234 56 : GEN Lk = gel(list,k), D = gel(Lk,3);
1235 56 : GEN sqrtD = sqrtr_abs(itor(D, prec)); /* sqrt(|D|) */
1236 56 : GEN indmultD = heegner_indexmultD(faN, indmult, itos(D), sqrtD);
1237 : do
1238 : {
1239 : GEN mulf, indr;
1240 : pari_timer ti;
1241 56 : if (DEBUGLEVEL) timer_start(&ti);
1242 56 : mulf = ltwist1(E, D, bprec+expo(indmultD));
1243 56 : if (DEBUGLEVEL) timer_printf(&ti,"ellL1twist");
1244 56 : indr = mulrr(indmultD, mulf);
1245 56 : if (DEBUGLEVEL) err_printf("Disc = %Ps, Index^2 = %Ps\n", D, indr);
1246 56 : if (signe(indr)>0 && expo(indr) >= -1) /* indr >=.5 */
1247 : {
1248 : long e, i, l;
1249 42 : GEN pts, cfs, L, indi = grndtoi(sqrtr_abs(indr), &e);
1250 42 : if (e > expi(indi)-7)
1251 : {
1252 0 : bprec++;
1253 0 : pari_warn(warnprec, "ellL1",bprec);
1254 0 : continue;
1255 : }
1256 42 : *pind = itos(indi);
1257 42 : L = gel(Lk, 2); l = lg(L);
1258 42 : pts = cgetg(l, t_VEC);
1259 42 : cfs = cgetg(l, t_VECSMALL);
1260 189 : for (i = 1; i < l; ++i)
1261 : {
1262 147 : GEN P = gel(L,i), z = gel(P,2), Q = gel(P,3); /* [1 or 2, Q] */
1263 : long c;
1264 147 : gel(pts, i) = qfb_root(gel(P,1), sqrtD);
1265 147 : c = rootno(Q, gel(faN,1), listR);
1266 147 : if (!equali1(z)) c *= 2;
1267 147 : cfs[i] = c;
1268 : }
1269 42 : if (DEBUGLEVEL)
1270 0 : err_printf("N = %Ps, ymin*N = %Ps\n",N,
1271 0 : gmul(gsqrt(gel(Lk, 1), prec),N));
1272 42 : *coefs = cfs; *points = pts; return;
1273 : }
1274 : } while(0);
1275 : }
1276 0 : d = listD[l-1]; set_avma(av);
1277 : }
1278 : }
1279 :
1280 : GEN
1281 154 : ellanal_globalred_all(GEN e, GEN *cb, GEN *N, GEN *tam)
1282 : {
1283 154 : GEN E = ellanal_globalred(e, cb), red = obj_check(E, Q_GLOBALRED);
1284 154 : *N = gel(red, 1);
1285 154 : *tam = gel(red,2);
1286 154 : if (signe(ell_get_disc(E))>0) *tam = shifti(*tam,1);
1287 154 : return E;
1288 : }
1289 :
1290 : static GEN
1291 42 : vecelnfembed(GEN x, GEN M, GEN et)
1292 91 : { pari_APPLY_same(gmul(M, etnf_to_basis(et, gel(x,i)))) }
1293 :
1294 : static GEN
1295 42 : QXQV_inv(GEN x, GEN T)
1296 91 : { pari_APPLY_same(QXQ_inv(gel(x,i), T)) }
1297 :
1298 : static GEN
1299 42 : etnfnewprec(GEN x, long prec)
1300 126 : { pari_APPLY_same(nfnewprec(gel(x,i),prec)) }
1301 :
1302 : static GEN
1303 49 : vec_etnf_to_basis(GEN et, GEN x)
1304 154 : { pari_APPLY_same(etnf_to_basis(et,gel(x,i))) }
1305 :
1306 : static GEN
1307 42 : makenfA(GEN sel, GEN A, GEN cb)
1308 : {
1309 42 : GEN etal = gel(sel,1), T = gel(etal,3);
1310 42 : GEN et = gel(etal,1), M = etnf_get_M(et);
1311 42 : long v = etnf_get_varn(et);
1312 42 : GEN al = vecelnfembed(A, M, et);
1313 42 : GEN th = gmul(M, etnf_to_basis(et, pol_x(v)));
1314 42 : return mkvec5(etal,QXQV_inv(A, T),cb,al,th);
1315 : }
1316 :
1317 : GEN
1318 56 : ellheegner(GEN E)
1319 : {
1320 56 : pari_sp av = avma;
1321 : GEN z, P, ht, points, coefs, s, om, indmult;
1322 : GEN sel, etal, et, cbb, A, dAi, T, Ag, At;
1323 : long ind, indx, lint, k, l, wtor, etor, ndisc, ltors2, selrank;
1324 56 : long bitprec = 16, prec = nbits2prec(bitprec) + EXTRAPRECWORD;
1325 : pari_timer ti;
1326 : GEN N, cb, tam, torsion, nfA;
1327 56 : E = ellanal_globalred_all(E, &cb, &N, &tam);
1328 56 : if (ellrootno_global(E) == 1)
1329 7 : pari_err_DOMAIN("ellheegner", "(analytic rank)%2","=",gen_0,E);
1330 49 : torsion = elltors(E);
1331 49 : wtor = itos( gel(torsion,1) ); /* #E(Q)_tor */
1332 49 : etor = wtor > 1? itou(gmael(torsion, 2, 1)): 1; /* exponent of E(Q)_tor */
1333 49 : sel = ell2selmer_basis(E, &cbb, prec);
1334 49 : etal = gel(sel,1); A = gel(sel,2); et = gel(etal,1); T = gel(etal,3);
1335 49 : ltors2 = lg(et)-2; selrank = lg(A)-1;
1336 49 : Ag = selrank > ltors2+1 ? pol_1(etnf_get_varn(et)): gel(A,selrank);
1337 49 : At = vecslice(A,1,ltors2);
1338 49 : dAi = gsupnorm(vec_etnf_to_basis(et,A),prec);
1339 : while (1)
1340 42 : {
1341 : GEN hnaive, l1;
1342 : long bitneeded;
1343 91 : if (DEBUGLEVEL) timer_start(&ti);
1344 91 : l1 = ellL1(E, 1, bitprec);
1345 91 : if (DEBUGLEVEL) timer_printf(&ti,"ellL1");
1346 91 : if (expo(l1) < 1 - bitprec/2)
1347 7 : pari_err_DOMAIN("ellheegner", "analytic rank",">",gen_1,E);
1348 84 : om = ellR_omega(E,prec);
1349 84 : ht = divrr(mulru(l1, wtor * wtor), mulri(gel(om,1), tam));
1350 84 : if (DEBUGLEVEL) err_printf("Expected height=%Ps\n", ht);
1351 84 : hnaive = hnaive_max(E, ht);
1352 84 : if (DEBUGLEVEL) err_printf("Naive height <= %Ps\n", hnaive);
1353 84 : hnaive = gadd(shiftr(hnaive,-1),glog(dAi,prec));
1354 84 : bitneeded = itos(gceil(divrr(hnaive, mplog2(prec)))) + 32;
1355 84 : if (DEBUGLEVEL) err_printf("precision = %ld\n", bitneeded);
1356 84 : if (bitprec>=bitneeded) break;
1357 42 : bitprec = bitneeded;
1358 42 : prec = nbits2prec(bitprec) + EXTRAPRECWORD;
1359 : }
1360 42 : indmult = heegner_indexmult(om, wtor, tam, prec);
1361 42 : ndisc = maxss(10, (long) rtodbl(ht)/10);
1362 42 : heegner_find_disc(&points, &coefs, &ind, E, indmult, ndisc, prec);
1363 42 : if (DEBUGLEVEL) timer_start(&ti);
1364 42 : s = heegner_psi(E, N, points, bitprec);
1365 42 : if (DEBUGLEVEL) timer_printf(&ti,"heegner_psi");
1366 42 : l = lg(points);
1367 42 : z = mulsr(coefs[1], gel(s, 1));
1368 147 : for (k = 2; k < l; ++k) z = addrr(z, mulsr(coefs[k], gel(s, k)));
1369 42 : z = gsub(z, gmul(gel(om,1), ground(gdiv(z, gel(om,1)))));
1370 42 : if (DEBUGLEVEL) err_printf("z=%.*Pg\n",nbits2ndec(bitprec), z);
1371 42 : lint = wtor > 1 ? ugcd(ind, etor): 1;
1372 42 : indx = lint*2*ind;
1373 42 : if (vals(indx) >= vals(etor))
1374 35 : A = mkvec(Ag);
1375 : else
1376 7 : A = mkvec2(Ag, QXQ_mul(Ag, gel(At,1), T));
1377 42 : gmael(sel,1,1) = etnfnewprec(et, prec);
1378 42 : nfA = makenfA(sel, A, cbb);
1379 42 : P = heegner_find_point(E, nfA, om, ht, gmulsg(2*lint, z), indx, prec);
1380 42 : if (DEBUGLEVEL) timer_printf(&ti,"heegner_find_point");
1381 42 : if (cb) P = ellchangepointinv(P, cb);
1382 42 : return gc_GEN(av, P);
1383 : }
1384 :
1385 : /* Modular degree */
1386 :
1387 : static GEN
1388 70 : ellisobound(GEN e)
1389 : {
1390 70 : GEN M = gel(ellisomat(e,0,1),2);
1391 70 : return vecmax(gel(M,1));
1392 : }
1393 : /* 4Pi^2 / E.area */
1394 : static GEN
1395 140 : getA(GEN E, long prec) { return mpdiv(sqrr(Pi2n(1,prec)), ellR_area(E, prec)); }
1396 :
1397 : /* Modular degree of elliptic curve e over Q, assuming Manin constant = 1
1398 : * (otherwise multiply by square of Manin constant). */
1399 : GEN
1400 70 : ellmoddegree(GEN E)
1401 : {
1402 70 : pari_sp av = avma;
1403 : GEN N, tam, mc2, d;
1404 : long b;
1405 70 : E = ellanal_globalred_all(E, NULL, &N, &tam);
1406 70 : mc2 = sqri(ellisobound(E));
1407 70 : b = expi(mulii(N, mc2)) + maxss(0, expo(getA(E, LOWDEFAULTPREC))) + 16;
1408 : for(;;)
1409 0 : {
1410 70 : long prec = nbits2prec(b), e, s;
1411 70 : GEN deg = mulri(mulrr(lfunellmfpeters(E, b), getA(E, prec)), mc2);
1412 70 : d = grndtoi(deg, &e);
1413 70 : if (DEBUGLEVEL) err_printf("ellmoddegree: %Ps, bit=%ld, err=%ld\n",deg,b,e);
1414 70 : s = expo(deg);
1415 70 : if (e <= -8 && s <= b-8) return gc_upto(av, gdiv(d,mc2));
1416 0 : b = maxss(s, b+e) + 16;
1417 : }
1418 : }
|