Testing latest pari + WASM + node.js... and it works?! Wow.
License: GPL3
ubuntu2004
/* Copyright (C) 2000 The PARI group.12This file is part of the PARI/GP package.34PARI/GP is free software; you can redistribute it and/or modify it under the5terms of the GNU General Public License as published by the Free Software6Foundation; either version 2 of the License, or (at your option) any later7version. It is distributed in the hope that it will be useful, but WITHOUT8ANY WARRANTY WHATSOEVER.910Check the License for details. You should have received a copy of it, along11with the package; see the file 'COPYING'. If not, write to the Free Software12Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */13#include "pari.h"14#include "paripriv.h"1516#define DEBUGLEVEL DEBUGLEVEL_hensel1718/***********************************************************************/19/** **/20/** QUADRATIC HENSEL LIFT (adapted from V. Shoup's NTL) **/21/** **/22/***********************************************************************/23/* Setup for divide/conquer quadratic Hensel lift24* a = set of k t_POL in Z[X] = factors over Fp (T=NULL) or Fp[Y]/(T)25* V = set of products of factors built as follows26* 1) V[1..k] = initial a27* 2) iterate:28* append to V the two smallest factors (minimal degree) in a, remove them29* from a and replace them by their product [net loss for a = 1 factor]30*31* W = bezout coeffs W[i]V[i] + W[i+1]V[i+1] = 132*33* link[i] = -j if V[i] = a[j]34* j if V[i] = V[j] * V[j+1]35* Arrays (link, V, W) pre-allocated for 2k - 2 elements */36static void37BuildTree(GEN link, GEN V, GEN W, GEN a, GEN T, GEN p)38{39long k = lg(a)-1;40long i, j, s, minp, mind;4142for (i=1; i<=k; i++) { gel(V,i) = gel(a,i); link[i] = -i; }4344for (j=1; j <= 2*k-5; j+=2,i++)45{46minp = j;47mind = degpol(gel(V,j));48for (s=j+1; s<i; s++)49if (degpol(gel(V,s)) < mind) { minp = s; mind = degpol(gel(V,s)); }5051swap(gel(V,j), gel(V,minp));52lswap(link[j], link[minp]);5354minp = j+1;55mind = degpol(gel(V,j+1));56for (s=j+2; s<i; s++)57if (degpol(gel(V,s)) < mind) { minp = s; mind = degpol(gel(V,s)); }5859swap(gel(V,j+1), gel(V,minp));60lswap(link[j+1], link[minp]);6162gel(V,i) = FqX_mul(gel(V,j), gel(V,j+1), T, p);63link[i] = j;64}6566for (j=1; j <= 2*k-3; j+=2)67{68GEN d, u, v;69d = FqX_extgcd(gel(V,j), gel(V,j+1), T, p, &u, &v);70if (degpol(d) > 0) pari_err_COPRIME("BuildTree", gel(V,j), gel(V,j+1));71d = gel(d,2);72if (!gequal1(d))73{74if (typ(d)==t_POL)75{76d = FpXQ_inv(d, T, p);77u = FqX_Fq_mul(u, d, T, p);78v = FqX_Fq_mul(v, d, T, p);79}80else81{82d = Fp_inv(d, p);83u = FqX_Fp_mul(u, d, T,p);84v = FqX_Fp_mul(v, d, T,p);85}86}87gel(W,j) = u;88gel(W,j+1) = v;89}90}9192/* au + bv = 1 (p0), ab = f (p0). Lift mod p1 = p0 pd (<= p0^2).93* If noinv is set, don't lift the inverses u and v */94static void95ZpX_HenselLift(GEN V, GEN W, long j, GEN f, GEN pd, GEN p0, GEN p1, int noinv)96{97pari_sp av = avma;98long space = lg(f) * lgefint(p1);99GEN a2, b2, g, z, s, t;100GEN a = gel(V,j), b = gel(V,j+1);101GEN u = gel(W,j), v = gel(W,j+1);102103(void)new_chunk(space); /* HACK */104g = ZX_sub(f, ZX_mul(a,b));105g = ZX_Z_divexact(g, p0);106g = FpX_red(g, pd);107z = FpX_mul(v,g, pd);108t = FpX_divrem(z,a, pd, &s);109t = ZX_add(ZX_mul(u,g), ZX_mul(t,b));110t = FpX_red(t, pd);111t = ZX_Z_mul(t,p0);112s = ZX_Z_mul(s,p0);113set_avma(av);114a2 = ZX_add(a,s);115b2 = ZX_add(b,t);116117/* already reduced mod p1 = pd p0 */118gel(V,j) = a2;119gel(V,j+1) = b2;120if (noinv) return;121122av = avma;123(void)new_chunk(space); /* HACK */124g = ZX_add(ZX_mul(u,a2), ZX_mul(v,b2));125g = Z_ZX_sub(gen_1, g);126g = ZX_Z_divexact(g, p0);127g = FpX_red(g, pd);128z = FpX_mul(v,g, pd);129t = FpX_divrem(z,a, pd, &s);130t = ZX_add(ZX_mul(u,g), ZX_mul(t,b));131t = FpX_red(t, pd);132t = ZX_Z_mul(t,p0);133s = ZX_Z_mul(s,p0);134set_avma(av);135gel(W,j) = ZX_add(u,t);136gel(W,j+1) = ZX_add(v,s);137}138139static void140ZpXQ_HenselLift(GEN V, GEN W, long j, GEN f, GEN Td, GEN T1, GEN pd, GEN p0, GEN p1, int noinv)141{142pari_sp av = avma;143const long n = degpol(T1), vT = varn(T1);144long space = lg(f) * lgefint(p1) * lg(T1);145GEN a2, b2, g, z, s, t;146GEN a = gel(V,j), b = gel(V,j+1);147GEN u = gel(W,j), v = gel(W,j+1);148149(void)new_chunk(space); /* HACK */150g = RgX_sub(f, Kronecker_to_ZXX(ZXX_mul_Kronecker(a,b,n), n, vT));151g = FpXQX_red(g, T1, p1);152g = RgX_Rg_divexact(g, p0);153z = FpXQX_mul(v,g, Td,pd);154t = FpXQX_divrem(z,a, Td,pd, &s);155t = ZX_add(ZXX_mul_Kronecker(u,g,n), ZXX_mul_Kronecker(t,b,n));156t = Kronecker_to_ZXX(t, n, vT);157t = FpXQX_red(t, Td, pd);158t = RgX_Rg_mul(t,p0);159s = RgX_Rg_mul(s,p0);160set_avma(av);161162a2 = RgX_add(a,s);163b2 = RgX_add(b,t);164/* already reduced mod p1 = pd p0 */165gel(V,j) = a2;166gel(V,j+1) = b2;167if (noinv) return;168169av = avma;170(void)new_chunk(space); /* HACK */171g = ZX_add(ZXX_mul_Kronecker(u,a2,n), ZXX_mul_Kronecker(v,b2,n));172g = Kronecker_to_ZXX(g, n, vT);173g = Rg_RgX_sub(gen_1, g);174g = FpXQX_red(g, T1, p1);175g = RgX_Rg_divexact(g, p0);176z = FpXQX_mul(v,g, Td,pd);177t = FpXQX_divrem(z,a, Td,pd, &s);178t = ZX_add(ZXX_mul_Kronecker(u,g,n), ZXX_mul_Kronecker(t,b,n));179t = Kronecker_to_ZXX(t, n, vT);180t = FpXQX_red(t, Td, pd);181t = RgX_Rg_mul(t,p0);182s = RgX_Rg_mul(s,p0);183set_avma(av);184gel(W,j) = RgX_add(u,t);185gel(W,j+1) = RgX_add(v,s);186}187188/* v list of factors, w list of inverses. f = v[j] v[j+1]189* Lift v[j] and v[j+1] mod p0 pd (possibly mod T), then all their divisors */190static void191ZpX_RecTreeLift(GEN link, GEN v, GEN w, GEN pd, GEN p0, GEN p1,192GEN f, long j, int noinv)193{194if (j < 0) return;195ZpX_HenselLift(v, w, j, f, pd, p0,p1, noinv);196ZpX_RecTreeLift(link, v, w, pd, p0,p1, gel(v,j) , link[j ], noinv);197ZpX_RecTreeLift(link, v, w, pd, p0,p1, gel(v,j+1), link[j+1], noinv);198}199static void200ZpXQ_RecTreeLift(GEN link, GEN v, GEN w, GEN Td, GEN T1, GEN pd, GEN p0, GEN p1,201GEN f, long j, int noinv)202{203if (j < 0) return;204ZpXQ_HenselLift(v, w, j, f, Td,T1, pd, p0,p1, noinv);205ZpXQ_RecTreeLift(link, v, w, Td,T1, pd, p0,p1, gel(v,j) , link[j ], noinv);206ZpXQ_RecTreeLift(link, v, w, Td,T1, pd, p0,p1, gel(v,j+1), link[j+1], noinv);207}208209/* Assume n > 0. We want to go to accuracy n, starting from accuracy 1, using210* a quadratically convergent algorithm. Goal: 9 -> 1,2,3,5,9 instead of211* 1,2,4,8,9 (sequence of accuracies).212*213* Let a0 = 1, a1 = 2, a2, ... ak = n, the sequence of accuracies. To obtain214* it, work backwards:215* a(k) = n, a(i-1) = (a(i) + 1) \ 2,216* but we do not want to store a(i) explicitly, even as a t_VECSMALL, since217* this would leave an object on the stack. We store a(i) implicitly in a218* MASK: let a(0) = 1, if the i-bit of MASK is set, set a(i+1) = 2 a(i) - 1,219* and 2a(i) otherwise.220*221* In fact, we do something a little more complicated to simplify the222* function interface and avoid returning k and MASK separately: we return223* MASK + 2^(k+1), so the highest bit of the mask indicates the length of the224* sequence, and the following ones are as above. */225ulong226quadratic_prec_mask(long n)227{228long a = n, i;229ulong mask = 0;230for(i = 1;; i++, mask <<= 1)231{232mask |= (a&1); a = (a+1)>>1;233if (a==1) return mask | (1UL << i);234}235}236237/* Lift to precision p^e0.238* a = modular factors of f mod (p,T) [possibly T=NULL]239* OR a TreeLift structure [e, link, v, w]: go on lifting240* flag = 0: standard.241* flag = 1: return TreeLift structure */242static GEN243MultiLift(GEN f, GEN a, GEN T, GEN p, long e0, long flag)244{245long i, eold, e, k = lg(a) - 1;246GEN E, v, w, link, penew, Tnew;247ulong mask;248pari_timer Ti;249250if (k < 2) pari_err_DOMAIN("MultiLift", "#(modular factors)", "<", gen_2,a);251if (e0 < 1) pari_err_DOMAIN("MultiLift", "precision", "<", gen_1,stoi(e0));252if (e0 == 1) return a;253254if (DEBUGLEVEL > 3) timer_start(&Ti);255if (typ(gel(a,1)) == t_INT)256{ /* a = TreeLift structure */257e = itos(gel(a,1));258link = gel(a,2);259v = gel(a,3);260w = gel(a,4);261}262else263{264e = 1;265v = cgetg(2*k-2 + 1, t_VEC);266w = cgetg(2*k-2 + 1, t_VEC);267link=cgetg(2*k-2 + 1, t_VECSMALL);268BuildTree(link, v, w, a, T? FpX_red(T,p): NULL, p);269if (DEBUGLEVEL > 3) timer_printf(&Ti, "building tree");270}271mask = quadratic_prec_mask(e0);272eold = 1;273penew = NULL;274Tnew = NULL;275if (DEBUGLEVEL > 3) err_printf("lifting to prec %ld\n", e0);276while (mask > 1)277{278long enew = eold << 1;279if (mask & 1) enew--;280mask >>= 1;281if (enew >= e) { /* mask == 1: last iteration */282GEN peold = penew? penew: powiu(p, eold);283GEN Td = NULL, pd;284long d = enew - eold; /* = eold or eold-1 */285/* lift from p^eold to p^enew */286pd = (d == eold)? peold: diviiexact(peold, p); /* p^d */287penew = mulii(peold,pd);288if (T) {289if (Tnew)290Td = (d == eold)? Tnew: FpX_red(Tnew,pd);291else292Td = FpX_red(T, peold);293Tnew = FpX_red(T, penew);294ZpXQ_RecTreeLift(link, v, w, Td, Tnew, pd, peold, penew, f, lg(v)-2,295(flag == 0 && mask == 1));296}297else298ZpX_RecTreeLift(link, v, w, pd, peold, penew, f, lg(v)-2,299(flag == 0 && mask == 1));300if (DEBUGLEVEL > 3) timer_printf(&Ti, "reaching prec %ld", enew);301}302eold = enew;303}304305if (flag)306E = mkvec4(utoipos(e0), link, v, w);307else308{309E = cgetg(k+1, t_VEC);310for (i = 1; i <= 2*k-2; i++)311{312long t = link[i];313if (t < 0) gel(E,-t) = gel(v,i);314}315}316return E;317}318319/* Q list of (coprime, monic) factors of pol mod (T,p). Lift mod p^e = pe.320* T may be NULL */321GEN322ZpX_liftfact(GEN pol, GEN Q, GEN pe, GEN p, long e)323{324pari_sp av = avma;325pol = FpX_normalize(pol, pe);326if (lg(Q) == 2) return mkvec(pol);327return gerepilecopy(av, MultiLift(pol, Q, NULL, p, e, 0));328}329330GEN331ZpXQX_liftfact(GEN pol, GEN Q, GEN T, GEN pe, GEN p, long e)332{333pari_sp av = avma;334pol = FpXQX_normalize(pol, T, pe);335if (lg(Q) == 2) return mkvec(pol);336return gerepilecopy(av, MultiLift(pol, Q, T, p, e, 0));337}338339GEN340ZqX_liftfact(GEN f, GEN a, GEN T, GEN pe, GEN p, long e)341{ return T ? ZpXQX_liftfact(f, a, T, pe, p, e): ZpX_liftfact(f, a, pe, p, e); }342GEN343ZqX_liftroot(GEN f, GEN a, GEN T, GEN p, long e)344{ return T ? ZpXQX_liftroot(f, a,T , p, e): ZpX_liftroot(f, a, p, e); }345346/* U = NULL treated as 1 */347static void348BezoutPropagate(GEN link, GEN v, GEN w, GEN pe, GEN U, GEN f, long j)349{350GEN Q, R;351if (j < 0) return;352353Q = FpX_mul(gel(v,j), gel(w,j), pe);354if (U)355{356Q = FpXQ_mul(Q, U, f, pe);357R = FpX_sub(U, Q, pe);358}359else360R = Fp_FpX_sub(gen_1, Q, pe);361gel(w,j+1) = Q; /* 0 mod U v[j], 1 mod (1-U) v[j+1] */362gel(w,j) = R; /* 1 mod U v[j], 0 mod (1-U) v[j+1] */363BezoutPropagate(link, v, w, pe, R, f, link[j ]);364BezoutPropagate(link, v, w, pe, Q, f, link[j+1]);365}366367/* as above, but return the Bezout coefficients for the lifted modular factors368* U[i] = 1 mod Qlift[i]369* 0 mod Qlift[j], j != i */370GEN371bezout_lift_fact(GEN pol, GEN Q, GEN p, long e)372{373pari_sp av = avma;374GEN E, link, v, w, pe;375long i, k = lg(Q)-1;376if (k == 1) retmkvec(pol_1(varn(pol)));377pe = powiu(p, e);378pol = FpX_normalize(pol, pe);379E = MultiLift(pol, Q, NULL, p, e, 1);380link = gel(E,2);381v = gel(E,3);382w = gel(E,4);383BezoutPropagate(link, v, w, pe, NULL, pol, lg(v)-2);384E = cgetg(k+1, t_VEC);385for (i = 1; i <= 2*k-2; i++)386{387long t = link[i];388if (t < 0) E[-t] = w[i];389}390return gerepilecopy(av, E);391}392393/* Front-end for ZpX_liftfact:394lift the factorization of pol mod p given by L to p^N (if possible) */395GEN396polhensellift(GEN pol, GEN L, GEN Tp, long N)397{398GEN T, p;399long i, l;400pari_sp av = avma;401void (*chk)(GEN, const char*);402403if (typ(pol) != t_POL) pari_err_TYPE("polhensellift",pol);404RgX_check_ZXX(pol, "polhensellift");405if (!is_vec_t(typ(L)) || lg(L) < 3) pari_err_TYPE("polhensellift",L);406if (N < 1) pari_err_DOMAIN("polhensellift", "precision", "<", gen_1,stoi(N));407if (!ff_parse_Tp(Tp, &T, &p, 0)) pari_err_TYPE("polhensellift",Tp);408chk = T? RgX_check_ZXX: RgX_check_ZX;409l = lg(L); L = leafcopy(L);410for (i = 1; i < l; i++)411{412GEN q = gel(L,i);413if (typ(q) != t_POL) gel(L,i) = scalar_ZX_shallow(q, varn(pol));414else chk(q, "polhensellift");415}416return gerepilecopy(av, ZqX_liftfact(pol, L, T, powiu(p,N), p, N));417}418419static GEN420FqV_roots_from_deg1(GEN x, GEN T, GEN p)421{422long i,l = lg(x);423GEN r = cgetg(l,t_COL);424for (i=1; i<l; i++) { GEN P = gel(x,i); gel(r,i) = Fq_neg(gel(P,2), T, p); }425return r;426}427428static GEN429ZpX_liftroots_full(GEN f, GEN S, GEN q, GEN p, long e)430{431pari_sp av = avma;432GEN y = ZpX_liftfact(f, deg1_from_roots(S, varn(f)), q, p, e);433return gerepileupto(av, FqV_roots_from_deg1(y, NULL, q));434}435436GEN437ZpX_roots(GEN F, GEN p, long e)438{439pari_sp av = avma;440GEN q = powiu(p, e);441GEN f = FpX_normalize(F, p);442GEN g = FpX_normalize(FpX_split_part(f, p), p);443GEN S;444if (degpol(g) < degpol(f))445{446GEN h = FpX_div(f, g, p);447F = gel(ZpX_liftfact(F, mkvec2(g, h), q, p, e), 1);448}449S = FpX_roots(g, p);450return gerepileupto(av, ZpX_liftroots_full(F, S, q, p, e));451}452453static GEN454ZpXQX_liftroots_full(GEN f, GEN S, GEN T, GEN q, GEN p, long e)455{456pari_sp av = avma;457GEN y = ZpXQX_liftfact(f, deg1_from_roots(S, varn(f)), T, q, p, e);458return gerepileupto(av, FqV_roots_from_deg1(y, T, q));459}460461GEN462ZpXQX_roots(GEN F, GEN T, GEN p, long e)463{464pari_sp av = avma;465GEN q = powiu(p, e);466GEN f = FpXQX_normalize(F, T, p);467GEN g = FpXQX_normalize(FpXQX_split_part(f, T, p), T, p);468GEN S;469if (degpol(g) < degpol(f))470{471GEN h = FpXQX_div(f, g, T, p);472F = gel(ZpXQX_liftfact(F, mkvec2(g, h), T, q, p, e), 1);473}474S = FpXQX_roots(g, T, p);475return gerepileupto(av, ZpXQX_liftroots_full(F, S, T, q, p, e));476}477478GEN479ZqX_roots(GEN F, GEN T, GEN p, long e)480{481return T ? ZpXQX_roots(F, T, p, e): ZpX_roots(F, p, e);482}483484/* SPEC:485* p is a t_INT > 1, e >= 1486* f is a ZX with leading term prime to p.487* a is a simple root mod l for all l|p.488* Return roots of f mod p^e, as integers (implicitly mod p^e)489* STANDARD USE: p is a prime power */490GEN491ZpX_liftroot(GEN f, GEN a, GEN p, long e)492{493pari_sp av = avma;494GEN q = p, fr, W;495ulong mask;496497a = modii(a,q);498if (e == 1) return a;499mask = quadratic_prec_mask(e);500fr = FpX_red(f,q);501W = Fp_inv(FpX_eval(ZX_deriv(fr), a, q), q); /* 1/f'(a) mod p */502for(;;)503{504q = sqri(q);505if (mask & 1) q = diviiexact(q, p);506mask >>= 1;507fr = FpX_red(f,q);508a = Fp_sub(a, Fp_mul(W, FpX_eval(fr, a,q), q), q);509if (mask == 1) return gerepileuptoint(av, a);510W = Fp_sub(shifti(W,1), Fp_mul(Fp_sqr(W,q), FpX_eval(ZX_deriv(fr),a,q), q), q);511}512}513514GEN515ZpX_liftroots(GEN f, GEN S, GEN p, long e)516{517long i, n = lg(S)-1, d = degpol(f);518GEN r;519if (n == d) return ZpX_liftroots_full(f, S, powiu(p, e), p, e);520r = cgetg(n+1, typ(S));521for (i=1; i <= n; i++)522gel(r,i) = ZpX_liftroot(f, gel(S,i), p, e);523return r;524}525526GEN527ZpXQX_liftroot_vald(GEN f, GEN a, long v, GEN T, GEN p, long e)528{529pari_sp av = avma, av2;530GEN pv = p, q, W, df, Tq, fr, dfr;531ulong mask;532a = Fq_red(a, T, p);533if (e <= v+1) return a;534df = RgX_deriv(f);535if (v) { pv = powiu(p,v); df = ZXX_Z_divexact(df, pv); }536mask = quadratic_prec_mask(e-v);537Tq = FpXT_red(T, p); dfr = FpXQX_red(df, Tq, p);538W = Fq_inv(FqX_eval(dfr, a, Tq, p), Tq, p); /* 1/f'(a) mod (T,p) */539q = p;540av2 = avma;541for (;;)542{543GEN u, fa, qv, q2v, q2, Tq2;544q2 = q; q = sqri(q);545if (mask & 1) q = diviiexact(q,p);546mask >>= 1;547if (v) { qv = mulii(q, pv); q2v = mulii(q2, pv); }548else { qv = q; q2v = q2; }549Tq2 = FpXT_red(T, q2v); Tq = FpXT_red(T, qv);550fr = FpXQX_red(f, Tq, qv);551fa = FqX_eval(fr, a, Tq, qv);552fa = typ(fa)==t_INT? diviiexact(fa,q2v): ZX_Z_divexact(fa, q2v);553a = Fq_sub(a, Fq_mul(Fq_mul(W,fa,Tq2,q2v),q2, Tq,qv), Tq, qv);554if (mask == 1) return gerepileupto(av, a);555dfr = FpXQX_red(df, Tq, q);556u = Fq_sub(Fq_mul(W,FqX_eval(dfr,a,Tq,q),Tq,q),gen_1,Tq,q);557u = typ(u)==t_INT? diviiexact(u,q2): ZX_Z_divexact(u,q2);558W = Fq_sub(W, Fq_mul(Fq_mul(u,W,Tq2,q2),q2, Tq,q),Tq,q);559if (gc_needed(av2,2))560{561if(DEBUGMEM>1) pari_warn(warnmem,"ZpXQX_liftroot, e = %ld", e);562gerepileall(av2, 3, &a, &W, &q);563}564}565}566567GEN568ZpXQX_liftroot(GEN f, GEN a, GEN T, GEN p, long e) { return ZpXQX_liftroot_vald(f,a,0,T,p,e); }569570GEN571ZpXQX_liftroots(GEN f, GEN S, GEN T, GEN p, long e)572{573long i, n = lg(S)-1, d = degpol(f);574GEN r;575if (n == d) return ZpXQX_liftroots_full(f, S, T, powiu(p, e), p, e);576r = cgetg(n+1, typ(S));577for (i=1; i <= n; i++)578gel(r,i) = ZpXQX_liftroot(f, gel(S,i), T, p, e);579return r;580}581582/* Same as ZpX_liftroot for the polynomial X^n-b*/583GEN584Zp_sqrtnlift(GEN b, GEN n, GEN a, GEN p, long e)585{586pari_sp ltop=avma;587GEN q, w, n_1;588ulong mask;589long pis2 = equalii(n, gen_2)? 1: 0;590if (e == 1) return icopy(a);591n_1 = subiu(n,1);592mask = quadratic_prec_mask(e);593w = Fp_inv(pis2 ? shifti(a,1): Fp_mul(n,Fp_pow(a,n_1,p), p), p);594q = p;595for(;;)596{597q = sqri(q);598if (mask & 1) q = diviiexact(q, p);599mask >>= 1;600if (lgefint(q) == 3 && lgefint(n) == 3)601{602ulong Q = uel(q,2), N = uel(n,2);603ulong A = umodiu(a, Q);604ulong B = umodiu(b, Q);605ulong W = umodiu(w, Q);606607A = Fl_sub(A, Fl_mul(W, Fl_sub(Fl_powu(A,N,Q), B, Q), Q), Q);608a = utoi(A);609if (mask == 1) break;610W = Fl_sub(Fl_add(W,W,Q),611Fl_mul(Fl_sqr(W,Q), Fl_mul(N,Fl_powu(A, N-1, Q), Q), Q), Q);612w = utoi(W);613}614else615{616/* a -= w (a^n - b) */617a = modii(subii(a, mulii(w, subii(Fp_pow(a,n,q),b))), q);618if (mask == 1) break;619/* w += w - w^2 n a^(n-1)*/620w = subii(shifti(w,1), Fp_mul(Fp_sqr(w,q),621pis2? shifti(a,1): mulii(n,Fp_pow(a,n_1,q)), q));622}623}624return gerepileuptoint(ltop,a);625}626627/* Same as ZpX_liftroot for the polynomial X^2-b */628GEN629Zp_sqrtlift(GEN b, GEN a, GEN p, long e)630{631return Zp_sqrtnlift(b, gen_2, a, p, e);632}633634GEN635Zp_sqrt(GEN x, GEN p, long e)636{637pari_sp av;638GEN z;639if (absequaliu(p,2)) return Z2_sqrt(x,e);640av = avma;641z = Fp_sqrt(Fp_red(x, p), p);642if (!z) return NULL;643if (e > 1) z = Zp_sqrtlift(x, z, p, e);644return gerepileuptoint(av, z);645}646647/* Compute (x-1)/(x+1)/p^k */648static GEN649ZpXQ_log_to_ath(GEN x, long k, GEN T, GEN p, long e, GEN pe)650{651pari_sp av = avma;652long vT = get_FpX_var(T);653GEN bn, bdi;654GEN bd = ZX_Z_add(x, gen_1);655if (absequaliu(p,2)) /*For p=2, we need to simplify by 2*/656{657bn = ZX_shifti(x,-(k+1));658bdi= ZpXQ_invlift(ZX_shifti(bd ,-1), pol_1(vT), T, p, e);659}660else661{662bn = ZX_Z_divexact(ZX_Z_sub(x, gen_1),powiu(p,k));663bdi= ZpXQ_invlift(bd, scalarpol(Fp_inv(gen_2,p),vT), T, p, e);664}665return gerepileupto(av, FpXQ_mul(bn, bdi, T, pe));666}667668/* Assume p odd, a = 1 [p], return log(a) */669GEN670ZpXQ_log(GEN a, GEN T, GEN p, long N)671{672pari_sp av = avma;673pari_timer ti;674long is2 = absequaliu(p,2);675ulong pp = is2 ? 0: itou_or_0(p);676double lp = is2 ? 1: pp ? log2(pp): expi(p);677long k = maxss(1 , (long) .5+pow((double)(N>>1)/(lp*lp), 1./3));678GEN ak, s, b, pol;679long e = is2 ? N-1: N;680long i, l = (e-2)/(2*(k+is2));681GEN pe = powiu(p,e);682GEN TNk, pNk = powiu(p,N+k);683if( DEBUGLEVEL>=3) timer_start(&ti);684TNk = FpX_get_red(get_FpX_mod(T), pNk);685ak = FpXQ_pow(a, powiu(p,k), TNk, pNk);686if( DEBUGLEVEL>=3) timer_printf(&ti,"FpXQ_pow(%ld)",k);687b = ZpXQ_log_to_ath(ak, k, T, p, e, pe);688if( DEBUGLEVEL>=3) timer_printf(&ti,"ZpXQ_log_to_ath");689pol= cgetg(l+3,t_POL);690pol[1] = evalsigne(1)|evalvarn(0);691for(i=0; i<=l; i++)692{693GEN g;694ulong z = 2*i+1;695if (pp)696{697long w = u_lvalrem(z, pp, &z);698g = powuu(pp,2*i*k-w);699}700else g = powiu(p,2*i*k);701gel(pol,i+2) = Fp_divu(g, z,pe);702}703if( DEBUGLEVEL>=3) timer_printf(&ti,"pol(%ld)",l);704s = FpX_FpXQ_eval(pol, FpXQ_sqr(b, T, pe), T, pe);705if( DEBUGLEVEL>=3) timer_printf(&ti,"FpX_FpXQ_eval");706s = ZX_shifti(FpXQ_mul(b, s, T, pe), 1);707return gerepileupto(av, is2? s: FpX_red(s, pe));708}709710/***********************************************************************/711/** **/712/** Generic quadratic hensel lift over Zp[X] **/713/** **/714/***********************************************************************/715/* q = p^N */716717GEN718gen_ZpM_Dixon(GEN F, GEN V, GEN q, GEN p, long N, void *E,719GEN lin(void *E, GEN F, GEN d, GEN q),720GEN invl(void *E, GEN d))721{722pari_sp av = avma;723long N2, M;724GEN VN2, V2, VM, bil;725GEN q2, qM;726V = FpM_red(V, q);727if (N == 1) return invl(E, V);728N2 = (N + 1)>>1; M = N - N2;729F = FpM_red(F, q);730qM = powiu(p, M);731q2 = M == N2? qM: mulii(qM, p);732/* q2 = p^N2, qM = p^M, q = q2 * qM */733VN2 = gen_ZpM_Dixon(F, V, q2, p, N2, E, lin, invl);734bil = lin(E, F, VN2, q);735V2 = ZM_Z_divexact(ZM_sub(V, bil), q2);736VM = gen_ZpM_Dixon(F, V2, qM, p, M, E, lin, invl);737return gerepileupto(av, FpM_red(ZM_add(VN2, ZM_Z_mul(VM, q2)), q));738}739740GEN741gen_ZpX_Dixon(GEN F, GEN V, GEN q, GEN p, long N, void *E,742GEN lin(void *E, GEN F, GEN d, GEN q),743GEN invl(void *E, GEN d))744{745pari_sp av = avma;746long N2, M;747GEN VN2, V2, VM, bil;748GEN q2, qM;749V = FpX_red(V, q);750if (N == 1) return invl(E, V);751N2 = (N + 1)>>1; M = N - N2;752F = FpXT_red(F, q);753qM = powiu(p, M);754q2 = M == N2? qM: mulii(qM, p);755/* q2 = p^N2, qM = p^M, q = q2 * qM */756VN2 = gen_ZpX_Dixon(F, V, q2, p, N2, E, lin, invl);757bil = lin(E, F, VN2, q);758V2 = ZX_Z_divexact(ZX_sub(V, bil), q2);759VM = gen_ZpX_Dixon(F, V2, qM, p, M, E, lin, invl);760return gerepileupto(av, FpX_red(ZX_add(VN2, ZX_Z_mul(VM, q2)), q));761}762763GEN764gen_ZpM_Newton(GEN x, GEN p, long n, void *E,765GEN eval(void *E, GEN f, GEN q),766GEN invd(void *E, GEN V, GEN v, GEN q, long M))767{768pari_sp ltop = avma, av;769long N = 1, N2, M;770long mask;771GEN q = p;772if (n == 1) return gcopy(x);773mask = quadratic_prec_mask(n);774av = avma;775while (mask > 1)776{777GEN qM, q2, v, V;778N2 = N; N <<= 1;779q2 = q;780if (mask&1UL) { /* can never happen when q2 = p */781N--; M = N2-1;782qM = diviiexact(q2,p); /* > 1 */783q = mulii(qM,q2);784} else {785M = N2;786qM = q2;787q = sqri(q2);788}789/* q2 = p^N2, qM = p^M, q = p^N = q2 * qM */790mask >>= 1;791v = eval(E, x, q);792V = ZM_Z_divexact(gel(v,1), q2);793x = FpM_sub(x, ZM_Z_mul(invd(E, V, v, qM, M), q2), q);794if (gc_needed(av, 1))795{796if(DEBUGMEM>1) pari_warn(warnmem,"gen_ZpX_Newton");797gerepileall(av, 2, &x, &q);798}799}800return gerepileupto(ltop, x);801}802803static GEN804_ZpM_invd(void *E, GEN V, GEN v, GEN q, long M/*unused*/)805{806(void)E; (void)M;807return FpM_mul(V, gel(v,2), q);808}809810static GEN811_ZpM_eval(void *E, GEN x, GEN q)812{813GEN f = RgM_Rg_sub_shallow(FpM_mul(x, FpM_red((GEN) E, q), q), gen_1);814return mkvec2(f, x);815}816817GEN818ZpM_invlift(GEN M, GEN C, GEN p, long n)819{820return gen_ZpM_Newton(C, p, n, (void *)M, _ZpM_eval, _ZpM_invd);821}822823GEN824gen_ZpX_Newton(GEN x, GEN p, long n, void *E,825GEN eval(void *E, GEN f, GEN q),826GEN invd(void *E, GEN V, GEN v, GEN q, long M))827{828pari_sp ltop = avma, av;829long N = 1, N2, M;830long mask;831GEN q = p;832if (n == 1) return gcopy(x);833mask = quadratic_prec_mask(n);834av = avma;835while (mask > 1)836{837GEN qM, q2, v, V;838N2 = N; N <<= 1;839q2 = q;840if (mask&1UL) { /* can never happen when q2 = p */841N--; M = N2-1;842qM = diviiexact(q2,p); /* > 1 */843q = mulii(qM,q2);844} else {845M = N2;846qM = q2;847q = sqri(q2);848}849/* q2 = p^N2, qM = p^M, q = p^N = q2 * qM */850mask >>= 1;851v = eval(E, x, q);852V = ZX_Z_divexact(gel(v,1), q2);853x = FpX_sub(x, ZX_Z_mul(invd(E, V, v, qM, M), q2), q);854if (gc_needed(av, 1))855{856if(DEBUGMEM>1) pari_warn(warnmem,"gen_ZpX_Newton");857gerepileall(av, 2, &x, &q);858}859}860return gerepileupto(ltop, x);861}862863struct _ZpXQ_inv864{865GEN T, a, p ,n;866};867868static GEN869_inv_invd(void *E, GEN V, GEN v, GEN q, long M/*unused*/)870{871struct _ZpXQ_inv *d = (struct _ZpXQ_inv *) E;872GEN Tq = FpXT_red(d->T, q);873(void)M;874return FpXQ_mul(V, gel(v,2), Tq, q);875}876877static GEN878_inv_eval(void *E, GEN x, GEN q)879{880struct _ZpXQ_inv *d = (struct _ZpXQ_inv *) E;881GEN Tq = FpXT_red(d->T, q);882GEN f = FpX_Fp_sub(FpXQ_mul(x, FpX_red(d->a, q), Tq, q), gen_1, q);883return mkvec2(f, x);884}885886GEN887ZpXQ_invlift(GEN a, GEN x, GEN T, GEN p, long e)888{889struct _ZpXQ_inv d;890d.a = a; d.T = T; d.p = p;891return gen_ZpX_Newton(x, p, e, &d, _inv_eval, _inv_invd);892}893894GEN895ZpXQ_inv(GEN a, GEN T, GEN p, long e)896{897pari_sp av=avma;898GEN ai;899if (lgefint(p)==3)900{901ulong pp = p[2];902ai = Flx_to_ZX(Flxq_inv(ZX_to_Flx(a,pp), ZXT_to_FlxT(T, pp), pp));903} else904ai = FpXQ_inv(FpX_red(a,p), FpXT_red(T,p),p);905return gerepileupto(av, ZpXQ_invlift(a, ai, T, p, e));906}907908GEN909ZpXQ_div(GEN a, GEN b, GEN T, GEN q, GEN p, long e)910{911return FpXQ_mul(a, ZpXQ_inv(b, T, p, e), T, q);912}913914GEN915ZpXQX_divrem(GEN x, GEN Sp, GEN T, GEN q, GEN p, long e, GEN *pr)916{917pari_sp av = avma;918GEN S = get_FpXQX_mod(Sp);919GEN b = leading_coeff(S), bi;920GEN S2, Q;921if (typ(b)==t_INT) return FpXQX_divrem(x, Sp, T, q, pr);922bi = ZpXQ_inv(b, T, p, e);923S2 = FqX_Fq_mul_to_monic(S, bi, T, q);924Q = FpXQX_divrem(x, S2, T, q, pr);925if (pr==ONLY_DIVIDES && !Q) { set_avma(av); return NULL; }926if (pr==ONLY_REM || pr==ONLY_DIVIDES) return gerepileupto(av, Q);927Q = FpXQX_FpXQ_mul(Q, bi, T, q);928gerepileall(av, 2, &Q, pr);929return Q;930}931932GEN933ZpXQX_digits(GEN x, GEN B, GEN T, GEN q, GEN p, long e)934{935pari_sp av = avma;936GEN b = leading_coeff(B), bi;937GEN B2, P, V, W;938long i, lV;939if (typ(b)==t_INT) return FpXQX_digits(x, B, T, q);940bi = ZpXQ_inv(b, T, p, e);941B2 = FqX_Fq_mul_to_monic(B, bi, T, q);942V = FpXQX_digits(x, B2, T, q);943lV = lg(V)-1;944P = FpXQ_powers(bi, lV-1, T, q);945W = cgetg(lV+1, t_VEC);946for(i=1; i<=lV; i++)947gel(W, i) = FpXQX_FpXQ_mul(gel(V,i), gel(P, i), T, q);948return gerepileupto(av, W);949}950951struct _ZpXQ_sqrtn952{953GEN T, a, n, ai;954};955956static GEN957_sqrtn_invd(void *E, GEN V, GEN v, GEN q, long M)958{959struct _ZpXQ_sqrtn *d = (struct _ZpXQ_sqrtn *) E;960GEN Tq = FpX_red(d->T, q), aiq = FpX_red(d->ai, q);961(void)M;962return FpXQ_mul(FpXQ_mul(V, gel(v,2), Tq, q), aiq, Tq, q);963}964965static GEN966_sqrtn_eval(void *E, GEN x, GEN q)967{968struct _ZpXQ_sqrtn *d = (struct _ZpXQ_sqrtn *) E;969GEN Tq = FpX_red(d->T, q);970GEN f = FpX_sub(FpXQ_pow(x, d->n, Tq, q), d->a, q);971return mkvec2(f, x);972}973974GEN975ZpXQ_sqrtnlift(GEN a, GEN n, GEN x, GEN T, GEN p, long e)976{977struct _ZpXQ_sqrtn d;978d.a = a; d.T = T; d.n = n;979d.ai = ZpXQ_inv(ZX_Z_mul(a, n),T,p,(e+1)>>1);980return gen_ZpX_Newton(x, p, e, &d, _sqrtn_eval, _sqrtn_invd);981}982983static GEN984to_ZX(GEN a, long v) { return typ(a)==t_INT? scalarpol_shallow(a,v): a; }985986GEN987Zq_sqrtnlift(GEN a, GEN n, GEN x, GEN T, GEN p, long e)988{989return T? ZpXQ_sqrtnlift(to_ZX(a,varn(T)), n, to_ZX(x,varn(T)), T, p, e)990: Zp_sqrtnlift(a, n, x, p, e);991}992993GEN994ZpXQ_sqrt(GEN a, GEN T, GEN p, long e)995{996pari_sp av = avma;997GEN z = FpXQ_sqrt(FpX_red(a, p), T, p);998if (!z) return NULL;999if (e <= 1) return gerepileupto(av, z);1000return gerepileupto(av, ZpXQ_sqrtnlift(a, gen_2, z, T, p, e));1001}10021003GEN1004ZpX_ZpXQ_liftroot_ea(GEN P, GEN S, GEN T, GEN p, long n, void *E,1005GEN early(void *E, GEN x, GEN q))1006{1007pari_sp ltop = avma, av;1008long N, r;1009long mask;1010GEN q2, q, W, Q, Tq2, Tq, Pq;1011pari_timer ti;1012T = FpX_get_red(T, powiu(p, n));1013if (n == 1) return gcopy(S);1014mask = quadratic_prec_mask(n);1015av = avma;1016q2 = p; q = sqri(p); mask >>= 1; N = 2;1017if (DEBUGLEVEL > 3) timer_start(&ti);1018Tq = FpXT_red(T,q);1019Tq2 = FpXT_red(Tq,q2);1020Pq = FpX_red(P,q);1021W = FpXQ_inv(FpX_FpXQ_eval(FpX_deriv(P,q2), S, Tq2, q2), Tq2, q2);1022Q = ZX_Z_divexact(FpX_FpXQ_eval(Pq, S, Tq, q), q2);1023r = brent_kung_optpow(degpol(P), 4, 3);1024if (DEBUGLEVEL > 3)1025err_printf("ZpX_ZpXQ_liftroot: lifting to prec %ld\n",n);1026for (;;)1027{1028GEN H, Sq, Wq, Spow, dP, qq, Pqq, Tqq;1029H = FpXQ_mul(W, Q, Tq2, q2);1030Sq = FpX_sub(S, ZX_Z_mul(H, q2), q);1031if (DEBUGLEVEL > 3)1032timer_printf(&ti,"ZpX_ZpXQ_liftroot: reaching prec %ld",N);1033if (mask==1)1034return gerepileupto(ltop, Sq);1035if (early)1036{1037GEN Se = early(E, Sq, q);1038if (Se) return gerepileupto(ltop, Se);1039}1040qq = sqri(q); N <<= 1;1041if (mask&1UL) { qq = diviiexact(qq, p); N--; }1042mask >>= 1;1043Pqq = FpX_red(P, qq);1044Tqq = FpXT_red(T, qq);1045Spow = FpXQ_powers(Sq, r, Tqq, qq);1046Q = ZX_Z_divexact(FpX_FpXQV_eval(Pqq, Spow, Tqq, qq), q);1047dP = FpX_FpXQV_eval(FpX_deriv(Pq, q), FpXV_red(Spow, q), Tq, q);1048Wq = ZX_Z_divexact(FpX_Fp_sub(FpXQ_mul(W, dP, Tq, q), gen_1, q), q2);1049Wq = ZX_Z_mul(FpXQ_mul(W, Wq, Tq2, q2), q2);1050Wq = FpX_sub(W, Wq, q);1051S = Sq; W = Wq; q2 = q; q = qq; Tq2 = Tq; Tq = Tqq; Pq = Pqq;1052if (gc_needed(av, 1))1053{1054if(DEBUGMEM>1) pari_warn(warnmem,"ZpX_ZpXQ_Newton");1055gerepileall(av, 8, &S, &W, &Q, &Tq2, &Tq, &Pq, &q, &q2);1056}1057}1058}10591060GEN1061ZpX_ZpXQ_liftroot(GEN P, GEN S, GEN T, GEN p, long n)1062{1063return ZpX_ZpXQ_liftroot_ea(P, S, T, p, n, NULL, NULL);1064}10651066GEN1067ZpX_Frobenius(GEN T, GEN p, long e)1068{1069return ZpX_ZpXQ_liftroot(get_FpX_mod(T), FpX_Frobenius(T, p), T, p, e);1070}10711072GEN1073ZpXQM_prodFrobenius(GEN M, GEN T, GEN p, long e)1074{1075pari_sp av = avma;1076GEN xp = ZpX_Frobenius(T, p, e);1077GEN z = FpXQM_autsum(mkvec2(xp, M), get_FpX_degree(T), T, powiu(p,e));1078return gerepilecopy(av, gel(z,2));1079}10801081/* Canonical lift of polynomial */10821083static GEN _can_invl(void *E, GEN V) {(void) E; return V; }10841085static GEN _can_lin(void *E, GEN F, GEN V, GEN q)1086{1087GEN v = RgX_splitting(V, 3);1088(void) E;1089return FpX_sub(V,ZXV_dotproduct(v, F), q);1090}10911092static GEN1093_can_iter(void *E, GEN f, GEN q)1094{1095GEN h = RgX_splitting(f,3);1096GEN h1s = ZX_sqr(gel(h,1)), h2s = ZX_sqr(gel(h,2)), h3s = ZX_sqr(gel(h,3));1097GEN h12 = ZX_mul(gel(h,1), gel(h,2));1098GEN h13 = ZX_mul(gel(h,1), gel(h,3));1099GEN h23 = ZX_mul(gel(h,2), gel(h,3));1100GEN h1c = ZX_mul(gel(h,1), h1s);1101GEN h3c = ZX_mul(gel(h,3), h3s);1102GEN th = ZX_mul(ZX_sub(h2s,ZX_mulu(h13,3)),gel(h,2));1103GEN y = FpX_sub(f,ZX_add(RgX_shift_shallow(h3c,2),ZX_add(RgX_shift_shallow(th,1),h1c)),q);1104(void) E;1105return mkvecn(7,y,h1s,h2s,h3s,h12,h13,h23);1106}11071108static GEN1109_can_invd(void *E, GEN V, GEN v, GEN qM, long M)1110{1111GEN h1s=gel(v,2), h2s=gel(v,3), h3s=gel(v,4);1112GEN h12=gel(v,5), h13=gel(v,6), h23=gel(v,7);1113GEN F = mkvec3(ZX_sub(h1s,RgX_shift_shallow(h23,1)),RgX_shift_shallow(ZX_sub(h2s,h13),1),1114ZX_sub(RgX_shift_shallow(h3s,2),RgX_shift_shallow(h12,1)));1115(void)E;1116return gen_ZpX_Dixon(ZXV_Z_mul(F, utoi(3)), V, qM, utoi(3), M, NULL,1117_can_lin, _can_invl);1118}11191120static GEN1121F3x_frobeniuslift(GEN P, long n)1122{ return gen_ZpX_Newton(Flx_to_ZX(P),utoi(3), n, NULL, _can_iter, _can_invd); }11231124static GEN _can5_invl(void *E, GEN V) {(void) E; return V; }11251126static GEN _can5_lin(void *E, GEN F, GEN V, GEN q)1127{1128ulong p = *(ulong*)E;1129GEN v = RgX_splitting(V, p);1130return FpX_sub(V,ZXV_dotproduct(v, F), q);1131}11321133/* P(X,t) -> P(X*t^n,t) mod (t^p-1) */1134static GEN1135_shift(GEN P, long n, ulong p, long v)1136{1137long i, l=lg(P);1138GEN r = cgetg(l,t_POL); r[1] = P[1];1139for(i=2;i<l;i++)1140{1141long s = n*(i-2)%p;1142GEN ci = gel(P,i);1143if (typ(ci)==t_INT)1144gel(r,i) = monomial(ci, s, v);1145else1146gel(r,i) = RgX_rotate_shallow(ci, s, p);1147}1148return FpXX_renormalize(r, l);1149}11501151struct _can_mul1152{1153GEN T, q;1154ulong p;1155};11561157static GEN1158_can5_mul(void *E, GEN A, GEN B)1159{1160struct _can_mul *d = (struct _can_mul *)E;1161GEN a = gel(A,1), b = gel(B,1);1162long n = itos(gel(A,2));1163GEN bn = _shift(b, n, d->p, get_FpX_var(d->T));1164GEN c = FpXQX_mul(a, bn, d->T, d->q);1165return mkvec2(c, addii(gel(A,2), gel(B,2)));1166}11671168static GEN1169_can5_sqr(void *E, GEN A)1170{1171return _can5_mul(E,A,A);1172}11731174static GEN1175_can5_iter(void *E, GEN f, GEN q)1176{1177pari_sp av = avma;1178struct _can_mul D;1179ulong p = *(ulong*)E;1180long i, vT = fetch_var();1181GEN N, P, d, V, fs;1182D.q = q; D.T = ZX_Z_sub(pol_xn(p,vT),gen_1);1183D.p = p;1184fs = mkvec2(_shift(f, 1, p, vT), gen_1);1185N = gel(gen_powu_i(fs,p-1,(void*)&D,_can5_sqr,_can5_mul),1);1186N = ZXX_evalx0(FpXQX_red(N,polcyclo(p,vT),q));1187P = FpX_mul(N,f,q);1188P = RgX_deflate(P, p);1189d = RgX_splitting(N, p);1190V = cgetg(p+1,t_VEC);1191gel(V,1) = ZX_mulu(gel(d,1), p);1192for(i=2; i<= (long)p; i++)1193gel(V,i) = ZX_mulu(RgX_shift_shallow(gel(d,p+2-i), 1), p);1194(void)delete_var(); return gerepilecopy(av, mkvec2(ZX_sub(f,P),V));1195}11961197static GEN1198_can5_invd(void *E, GEN H, GEN v, GEN qM, long M)1199{1200ulong p = *(long*)E;1201return gen_ZpX_Dixon(gel(v,2), H, qM, utoipos(p), M, E, _can5_lin, _can5_invl);1202}12031204GEN1205Flx_Teichmuller(GEN P, ulong p, long n)1206{1207return p==3 ? F3x_frobeniuslift(P,n):1208gen_ZpX_Newton(Flx_to_ZX(P),utoipos(p), n, &p, _can5_iter, _can5_invd);1209}12101211GEN1212polteichmuller(GEN P, ulong p, long n)1213{1214pari_sp av = avma;1215GEN q = NULL;1216if (typ(P)!=t_POL || !RgX_is_FpX(P,&q)) pari_err_TYPE("polteichmuller",P);1217if (q && !equaliu(q,p)) pari_err_MODULUS("polteichmuller",q,utoi(p));1218if (n <= 0)1219pari_err_DOMAIN("polteichmuller", "precision", "<=",gen_0,stoi(n));1220return gerepileupto(av, p==2 ? F2x_Teichmuller(RgX_to_F2x(P), n)1221: Flx_Teichmuller(RgX_to_Flx(P, p), p, n));1222}122312241225