/***********************************************************************/ /** **/ /** ARITHMETIC OPERATIONS ON POLYNOMIALS **/ /** (first part) **/ /** **/ /***********************************************************************/ /* $Id: polarit1.c,v 1.2 1999/09/23 17:50:56 karim Exp $ */ #include "pari.h" GEN get_mul_table(GEN x,GEN bas,GEN *T); GEN pol_to_monic(GEN pol, GEN *lead); GEN sort_factor(GEN y, int (*cmp)(GEN,GEN)); GEN bsrch(GEN p, GEN fa, long Ka, GEN eta, long Ma); GEN eleval(GEN f,GEN h,GEN a); GEN respm(GEN f1,GEN f2,GEN pm); GEN setup(GEN p,GEN f,GEN theta,GEN nut, long *La, long *Ma); GEN vstar(GEN p,GEN h); /* see splitgen() for how to use these two */ GEN setloop(GEN a) { a=icopy(a); new_chunk(2); /* dummy to get two cells of extra space */ return a; } /* assume a > 0 */ GEN incpos(GEN a) { long i,l=lgefint(a); for (i=l-1; i>1; i--) if (++a[i]) return a; i=l+1; a--; /* use extra cell */ a[0]=evaltyp(1) | evallg(i); a[1]=evalsigne(1) | evallgefint(i); return a; } GEN incloop(GEN a) { long i,l; switch(signe(a)) { case 0: a--; /* use extra cell */ a[0]=evaltyp(t_INT) | evallg(3); a[1]=evalsigne(1) | evallgefint(3); a[2]=1; return a; case -1: l=lgefint(a); for (i=l-1; i>1; i--) if (a[i]--) break; if (a[2] == 0) { a++; /* save one cell */ a[0] = evaltyp(t_INT) | evallg(2); a[1] = evalsigne(0) | evallgefint(2); } return a; default: return incpos(a); } } /*******************************************************************/ /* */ /* DIVISIBILITE */ /* Renvoie 1 si y divise x, 0 sinon . */ /* */ /*******************************************************************/ int gdivise(GEN x, GEN y) { long av=avma; x=gmod(x,y); avma=av; return gcmp0(x); } int poldivis(GEN x, GEN y, GEN *z) { long av = avma; GEN p1 = poldivres(x,y,ONLY_DIVIDES); if (p1) { *z = p1; return 1; } avma=av; return 0; } /*******************************************************************/ /* */ /* REDUCTION */ /* Do the transformation t_FRACN/t_RFRACN --> t_FRAC/t_RFRAC */ /* */ /*******************************************************************/ /* x[1] is scalar, non-zero */ static GEN gred_simple(GEN x) { GEN p1,p2,x2,x3; x2=content((GEN)x[2]); if (gcmp1(x2)) return gcopy(x); x3=gdiv((GEN)x[1],x2); p2=denom(x3); x2=gdiv((GEN)x[2],x2); p1=cgetg(3,t_RFRAC); p1[1]=(long)numer(x3); p1[2]=lmul(x2,p2); return p1; } GEN gred_rfrac(GEN x) { GEN y,p1,xx1,xx2,x3, x1 = (GEN)x[1], x2 = (GEN)x[2]; long tx,ty; if (gcmp0(x1)) return gcopy(x1); tx=typ(x1); ty=typ(x2); if (ty!=t_POL) { if (tx!=t_POL) return gcopy(x); if (gvar2(x2) > varn(x1)) return gdiv(x1,x2); err(talker,"incompatible variables in gred"); } if (tx!=t_POL) { if (varn(x2) < gvar2(x1)) return gred_simple(x); err(talker,"incompatible variables in gred"); } if (varn(x2) < varn(x1)) return gred_simple(x); if (varn(x2) > varn(x1)) return gdiv(x1,x2); /* now x1 and x2 are polynomials with the same variable */ xx1=content(x1); if (!gcmp1(xx1)) x1=gdiv(x1,xx1); xx2=content(x2); if (!gcmp1(xx2)) x2=gdiv(x2,xx2); x3=gdiv(xx1,xx2); y = poldivres(x1,x2,&p1); if (!signe(p1)) { cgiv(p1); return gmul(x3,y); } p1 = ggcd(x2,p1); if (!isscalar(p1)) { x1=gdeuc(x1,p1); x2=gdeuc(x2,p1); } xx1=numer(x3); xx2=denom(x3); p1=cgetg(3,t_RFRAC); p1[1]=lmul(x1,xx1); p1[2]=lmul(x2,xx2); return p1; } GEN gred(GEN x) { long tx=typ(x),av=avma; GEN y,p1,x1,x2; if (is_frac_t(tx)) { x1=(GEN)x[1]; x2=(GEN)x[2]; y = dvmdii(x1,x2,&p1); if (p1 == gzero) return y; /* gzero volontaire */ (void)new_chunk((lgefint(x1)+lgefint(x2))<<1); p1=mppgcd(x2,p1); if (is_pm1(p1)) { avma=av; return gcopy(x); } avma=av; y=cgetg(3,t_FRAC); y[1]=ldivii(x1,p1); y[2]=ldivii(x2,p1); return y; } if (is_rfrac_t(tx)) return gerepileupto(av, gred_rfrac(x)); return gcopy(x); } /*******************************************************************/ /* */ /* POLYNOMIAL EUCLIDEAN DIVISION */ /* */ /*******************************************************************/ GEN gdivexact(GEN x, GEN y); /* Polynomial division x / y: * if z = ONLY_REM return remainder, otherwise return quotient * if z != NULL set *z to remainder * *z is the last object on stack (and thus can be disposed of with cgiv * instead of gerepile) */ GEN poldivres(GEN x, GEN y, GEN *pr) { long ty=typ(y),tx,vx,vy,dx,dy,dz,i,j,avy,av,av1,sx,lrem; int remainder; GEN z,p1,rem,y_lead,mod; GEN (*f)(GEN,GEN); if (pr == ONLY_DIVIDES_EXACT) { f = gdivexact; pr = ONLY_DIVIDES; } else f = gdiv; if (is_scalar_t(ty)) { if (pr == ONLY_REM) return gzero; if (pr && pr != ONLY_DIVIDES) *pr=gzero; return f(x,y); } tx=typ(x); vy=gvar9(y); if (is_scalar_t(tx) || gvar9(x)>vy) { if (pr == ONLY_REM) return gcopy(x); if (pr == ONLY_DIVIDES) return gcmp0(x)? gzero: NULL; if (pr) *pr=gcopy(x); return gzero; } if (tx!=t_POL || ty!=t_POL) err(typeer,"euclidean division (poldivres)"); vx=varn(x); if (vx=0; dy--) { y_lead = (GEN)y[dy+2]; if (!gcmp0(y_lead)) break; } } if (!dy) /* y is constant */ { if (pr && pr != ONLY_DIVIDES) { if (pr == ONLY_REM) return zeropol(vx); *pr = zeropol(vx); } return f(x,(GEN)y[2]); } dx=lgef(x)-3; if (vx>vy || dx=dy; i--) { av1=avma; p1=(GEN)x[i]; for (j=i-dy+1; j<=i && j<=dz; j++) if (y[i-j]) p1 = gadd(p1, gmul((GEN)z[j],(GEN)y[i-j])); if (y_lead) p1 = f(p1,y_lead); if (!remainder) p1 = avma==av1? gcopy(p1): gerepileupto(av1,p1); z[i-dy] = (long)p1; } if (!pr) return gerepileupto(av,z-2); rem = (GEN)avma; av1 = (long)new_chunk(dx+3); for (sx=0; ; i--) { p1 = (GEN)x[i]; /* we always enter this loop at least once */ for (j=0; j<=i && j<=dz; j++) if (y[i-j]) p1 = gadd(p1, gmul((GEN)z[j],(GEN)y[i-j])); if (mod && avma==av1) p1 = gmod(p1,mod); if (!gcmp0(p1)) { sx = 1; break; } /* remainder is non-zero */ if (!isinexactreal(p1) && !isexactzero(p1)) break; if (!i) break; avma=av1; } if (pr == ONLY_DIVIDES) { if (sx) { avma=av; return NULL; } avma = (long)rem; return gerepileupto(av,z-2); } lrem=i+3; rem -= lrem; if (avma==av1) { avma = (long)rem; p1 = gcopy(p1); } else p1 = gerepileupto((long)rem,p1); rem[0]=evaltyp(t_POL) | evallg(lrem); rem[1]=evalsigne(1) | evalvarn(vx) | evallgef(lrem); rem += 2; rem[i]=(long)p1; for (i--; i>=0; i--) { av1=avma; p1 = (GEN)x[i]; for (j=0; j<=i && j<=dz; j++) if (y[i-j]) p1 = gadd(p1, gmul((GEN)z[j],(GEN)y[i-j])); if (mod && avma==av1) p1 = gmod(p1,mod); rem[i]=avma==av1? lcopy(p1):lpileupto(av1,p1); } rem -= 2; if (!sx) normalizepol_i(rem, lrem); if (remainder) return gerepileupto(av,rem); z -= 2; { GEN *gptr[2]; gptr[0]=&z; gptr[1]=&rem; gerepilemanysp(av,avy,gptr,2); *pr = rem; return z; } } /*******************************************************************/ /* */ /* ROOTS MODULO a prime p (no multiplicities) */ /* */ /*******************************************************************/ static GEN mod(GEN x, GEN y) { GEN z = cgetg(3,t_INTMOD); z[1]=(long)y; z[2]=(long)x; return z; } static long factmod_init(GEN *F, GEN pp, long *p) { GEN f = *F; long i,d; if (typ(f)!=t_POL || typ(pp)!=t_INT) err(typeer,"factmod"); if (expi(pp) > BITS_IN_LONG - 3) *p = 0; else { *p = itos(pp); if (*p < 2) err(talker,"not a prime in factmod"); } f = gmul(f, mod(gun,pp)); if (!signe(f)) err(zeropoler,"factmod"); f = lift_intern(f); d = lgef(f); for (i=2; i ss[2]); if (nbrac == 1) { avma=av; return cgetg(1,t_COL); } if (nbrac == d && p != ss[2]) { g = mpinvmod((GEN)f[3], pp); setsigne(g,-1); ss = modis(mulii(g, (GEN)f[2]), p); y[nbrac++]=ss[2]; } avma=av; g=cgetg(nbrac,t_COL); if (isonstack(pp)) pp = icopy(pp); for (i=1; i1) y[1] = (long)gmodulsg(0,p); return y; } y = cgetg(n+j,t_COL); if (j>1) { y[1] = zero; n++; } y[j] = (long)normalize_mod_p(b,p); if (la) y[j+lb] = (long)normalize_mod_p(a,p); pol = gadd(polx[varn(f)], gun); pol0 = (GEN)pol[2]; while (j<=n) { a=(GEN)y[j]; la=lgef(a)-3; if (la==1) y[j++] = lsubii(p, (GEN)a[2]); else if (la==2) { GEN d = subii(sqri((GEN)a[3]), shifti((GEN)a[2],2)); GEN e = mpsqrtmod(d,p), u = addis(q, 1); /* u = 1/2 */ y[j++] = lmodii(mulii(u,subii(e,(GEN)a[3])), p); y[j++] = lmodii(mulii(u,negi(addii(e,(GEN)a[3]))), p); } else for (pol0[2]=1; ; pol0[2]++) { b = Fp_pow_mod_pol(pol,q, a,p); b[2] = laddis((GEN)b[2], -1); b = Fp_pol_gcd(a,b, p); lb = lgef(b)-3; if (lb && lb3) { k++; if (p && !(k%p)) { k++; f2 = Fp_deuc(f2,g1,pp); } p1 = Fp_pol_gcd(f2,g1, pp); u = g1; g1 = p1; if (!gcmp1(p1)) { u = Fp_deuc( u,p1,pp); f2= Fp_deuc(f2,p1,pp); } if (lgef(u)>3) { /* here u is square-free (product of irred. of multiplicity e * k) */ pd=gun; v=polx[vf]; for (d=1; d <= (lgef(u)-3)>>1; d++) { pd=mulii(pd,pp); v=Fp_pow_mod_pol(v,pp, u,pp); g = Fp_pol_gcd(gadd(v, gneg(polx[vf])), u, pp); if (lgef(g)>3) { /* Ici g est produit de pol irred ayant tous le meme degre d; */ j=nbfact+(lgef(g)-3)/d; if (simple) for ( ; nbfact X) pour faire pgcd de g avec * w^(p^d-1)/2 jusqu'a casser. p = 2 will be treated separately. */ if (p) split(p,t+nbfact,d,pp,q,r); else splitgen(pp,t+nbfact,d,pp,q,r); for (; nbfact3) { t[nbfact] = simple? (GEN)(lgef(u)-3): u; ex[nbfact++]=e*k; } } } j = lgef(f2); if (j==3) break; e*=p; j=(j-3)/p+3; setlg(f,j); setlgef(f,j); for (i=2; i1 && !x[i]); if (i==1) setsigne(x,0); else { setsigne(x,1); setlgef(x,i+1); } } long split_berlekamp(GEN Q, GEN *t, GEN pp, GEN pps2) { GEN u = *t, p1, p2, vker,v,w,pol; long av,N = lgef(u)-3, d,i,j,kk,l1,l2,p, vu = varn(u); if (DEBUGLEVEL > 7) timer2(); p = is_bigint(pp)? 0: pp[2]; setlg(Q, N+1); setlg(Q[1], N+1); w = v = Fp_pow_mod_pol(polx[vu],pp,u,pp); for (j=2; j<=N; j++) { p1 = (GEN)Q[j]; setlg(p1, N+1); d=lgef(w)-1; p2 = w+1; for (i=1; i 7) msgtimer("frobenius"); vker = mat_to_vecpol(ker_mod_p(Q,pp), vu); if (DEBUGLEVEL > 7) msgtimer("kernel"); d = lg(vker)-1; if (p) for (i=1; i<=d; i++) { p1 = (GEN)vker[i]; for (j=2; j1) { av = avma; p2 = Fp_res(polt, p1, pp); if (lgef(p2) <= 3) { avma=av; continue; } p2 = Fp_pow_mod_pol(p2,pps2, p1,pp); /* set p2 = p2 - 1 */ p2[2] = laddis((GEN)p2[2], -1); p2 = Fp_pol_gcd(p1,p2, pp); l2=lgef(p2)-3; if (l2>0 && l2 7) msgtimer("new factor"); } else avma = av; } } } return d; } GEN factmod(GEN f, GEN pp) { long i,j,k,e,p,N,nbfact,av = avma,tetpil,d; GEN pps2,ex,y,f2,p1,g1,Q,u,v, *t; if (!(d = factmod_init(&f, pp, &p))) { avma=av; return trivfact(); } /* to hold factors and exponents */ t = (GEN*)cgetg(d+1,t_VEC); ex = new_chunk(d+1); e = nbfact = 1; pps2 = shifti(pp,-1); Q = cgetg(d+1,t_MAT); for (i=1; i<=d; i++) Q[i] = lgetg(d+1, t_COL); p1 = (GEN)Q[1]; for (i=1; i<=d; i++) p1[i] = zero; for(;;) { f2 = Fp_pol_gcd(f,derivpol(f), pp); g1 = gcmp1(f2)? f: Fp_deuc(f,f2,pp); k = 0; while (lgef(g1)>3) { k++; if (p && !(k%p)) { k++; f2 = Fp_deuc(f2,g1,pp); } p1 = Fp_pol_gcd(f2,g1, pp); u = g1; g1 = p1; if (!gcmp1(p1)) { u = Fp_deuc( u,p1,pp); f2= Fp_deuc(f2,p1,pp); } N = lgef(u)-3; if (N) { /* here u is square-free (product of irred. of multiplicity e * k) */ t[nbfact] = normalize_mod_p(u,pp); d = (N==1)? 1: split_berlekamp(Q, t+nbfact, pp, pps2); for (j=0; j3) { f=gdeuc(f,p1); fp=derivpol(f); } p=(GEN)a[2]; p1=poleval(f,a); v=ggval(p1,p); if (v <= 0) err(rootper2); fl2=egalii(p,gdeux); if (fl2) { if (v==1) err(rootper2); quatre=stoi(4); } v=ggval(poleval(fp,a),p); if (!v) /* simple zero */ { while (!gcmp0(p1)) { a = gsub(a,gdiv(p1,poleval(fp,a))); p1 = poleval(f,a); } tetpil=avma; pro=cgetg(2,t_VEC); pro[1]=lcopy(a); return gerepile(av,tetpil,pro); } n=lgef(f)-3; pro=cgetg(n+1,t_VEC); p1=poleval(f,gadd(a,gmul(fl2?quatre:p,polx[varn(f)]))); if (!gcmp0(p1)) p1=gdiv(p1,gpuigs(p,ggval(p1,p))); if (is_bigint(p)) err(impl,"apprgen for p>=2^31"); x = ggrandocp(p, valp(a) | precp(a)); if (fl2) { ps=4; x2=ggrandocp(p,2); p=quatre; } else { ps=itos(p); x2=ggrandocp(p,1); } for (j=0,i=0; i3) { f=gdeuc(f,p1); fp=derivpol(f); } fl2=egalii(p,gdeux); rac=(fl2 && r>=2)? rootmod(f,stoi(4)): rootmod(f,p); lx=lg(rac); p=gclone(p); if (r==1) { tetpil=avma; y=cgetg(lx,t_COL); for (i=1; i2, r>=2. On rend les n racines p-adiques en precision r si flall>0, * 1 seule si flall=0 */ GEN rootpadicfast(GEN f, GEN p, long r, long flall) { GEN rac,fp,y,z,p1,p2; long i,e,n; rac=rootmod(f,p); n=flall? lgef(f)-3: 1; fp=derivpol(f); y=cgetg(n+1,t_COL); p=gclone(p); p2=NULL; z=cgetg(5,t_PADIC); z[2]=(long)p; for (i=1; i<=n; i++) { p1=gmael(rac,i,2); if (signe(p1)) { if (!p2) p2=sqri(p); z[1] = evalvalp(0)|evalprecp(2); z[3] = (long)p2; } else { z[1] = evalvalp(2); z[3] = un; } z[4] = (long)p1; p1 = z; for(e=2;;) { p1 = gsub(p1, gdiv(poleval(f,p1),poleval(fp,p1))); if (e==r) break; e<<=1; if (e>r) e=r; p1 = gprec(p1,e); } y[i] = (long)p1; } return y; } static long getprec(GEN x, long prec, GEN *p) { long i,e; GEN p1; for (i = lgef(x)-1; i>1; i--) { p1=(GEN)x[i]; if (typ(p1)==t_PADIC) { e=valp(p1); if (signe(p1[4])) e += precp(p1); if (e3) { f=gdeuc(f,p1); fp=derivpol(f); } t=(GEN)a[1]; prec = getprec((GEN)a[2], BIGINT, &p); prec = getprec(t, prec, &p); if (prec==BIGINT) err(rootper1); p1=poleval(f,a); v=ggval(lift_intern(p1),p); if (v<=0) err(rootper2); fl2=egalii(p,gdeux); if (fl2) { if (v==1) err(rootper2); quatre=stoi(4); } v=ggval(lift_intern(poleval(fp,a)), p); if (!v) { while (!gcmp0(p1)) { a = gsub(a, gdiv(p1,poleval(fp,a))); p1 = poleval(f,a); } tetpil=avma; pro=cgetg(2,t_COL); pro[1]=lcopy(a); return gerepile(av,tetpil,pro); } n=lgef(f)-3; pro=cgetg(n+1,t_COL); j=0; p1=poleval(f,gadd(a,gmul(fl2?quatre:p,polx[varn(f)]))); if (!gcmp0(p1)) p1=gdiv(p1,gpuigs(p,ggval(p1,p))); if (is_bigint(p)) err(impl,"apprgen9 for p>=2^31"); x=gmodulcp(ggrandocp(p,prec), t); if (fl2) { ps_1=3; x2=ggrandocp(p,2); p=quatre; } else { ps_1=itos(p)-1; x2=ggrandocp(p,1); } d=lgef(t)-3; vecg=cgetg(d+1,t_COL); for (i=1; i<=d; i++) vecg[i] = (long)setloop(gzero); va=varn(t); for(;;) /* loop through F_q */ { ip=gmodulcp(gtopoly(vecg,va),t); if (gcmp0(poleval(p1,gadd(ip,x2)))) { u=apprgen9(p1,gadd(ip,x)); lu=lg(u); for (k=1; k3) ? padicff(res,p,r) : cgetg(1,t_COL); nbfac += (lg(fa[j])-1); } av2=avma; y=cgetg(3,t_MAT); p1=cgetg(nbfac+1,t_COL); y[1]=(long)p1; p2=cgetg(nbfac+1,t_COL); y[2]=(long)p2; for (i=1,k=0; i<=j; i++) for (i1=1; i1=3) { fprintferr(" entering Decomp_padic "); if (DEBUGLEVEL>=4) { fprintferr("with params: p=%Z, exponent=%ld, prec=%ld\n", p,mf,r); fprintferr(" f=%Z",f); } fprintferr("\n"); } unmodp=gmodulsg(1,p); pdr = respm(f, derivpol(f), gpuigs(p,mf)); b1=lift_intern(gmul(chi,unmodp)); a1=gun; b2=gun; b3=lift_intern(gmul(nu,unmodp)); while (lgef(b3) > 3) { GEN p1; b1 = Fp_deuc(b1,b3, p); b2 = Fp_pol_red(gmul(b2,b3), p); b3 = Fp_pol_extgcd(b2,b1, p, &a1,&p1); p1 = leading_term(b3); if (!gcmp1(p1)) { p1 = mpinvmod(p1,p); b3 = gmul(b3,p1); a1 = gmul(a1,p1); } } e=eleval(f,Fp_pol_red(gmul(a1,b2), p),theta); if (padicprec(e,p) > 0) e=gdiv(polmodi(gmul(pdr,e), mulii(pdr,p)),pdr); pk=p; pr=gpuigs(p,r); ph=mulii(pdr,pr); valk = 1; /* E(t)-e(t) belongs to p^k Op, which is contained in p^(k-df)*Zp[xi] */ while (cmpii(pk,ph) < 0) { e = gmul(e,gmul(e,gsubsg(3,gmul2n(e,1)))); e = gres(e,f); pk=sqri(pk); valk <<= 1; if (valk<=padicprec(e,p)) e = gdiv(polmodi(gmul(pdr,e), mulii(pk,pdr)), pdr); } f1=gcdpm(f,gmul(pdr,gsubsg(1,e)), ph); f1 = Fp_res(f1,f, pr); f2 = Fp_res(Fp_deuc(f,f1, pr), f, pr); if (DEBUGLEVEL>=4) { fprintferr(" leaving Decomp_padic with parameters: "); fprintferr("f1=%Z, f2=%Z\n",f1,f2); } b1=factorpadic4(f1,p,r); b2=factorpadic4(f2,p,r); res=cgetg(3,t_MAT); res[1]=lconcat((GEN)b1[1],(GEN)b2[1]); res[2]=lconcat((GEN)b1[2],(GEN)b2[2]); return res; } static GEN nilordpadic(GEN p,long r,GEN fx,long mf,GEN gx) { long Da,Na,La,Ma,first=1,n,v=varn(fx); GEN alpha,chi,nu,eta,w,phi,pmf,Dchi; if (DEBUGLEVEL>=3) { fprintferr(" entering Nilord_padic "); if (DEBUGLEVEL>=4) { fprintferr("with parameters: p=%Z, expo=%ld\n",p,mf); fprintferr(" fx=%Z, gx=%Z",fx,gx); } fprintferr("\n"); } pmf=gpuigs(p,mf+1); alpha=polx[v]; n=lgef(fx)-3; for(;;) { if (first) { chi=fx; nu=gx; first=0; } else { w=factcp(p,fx,alpha); chi=(GEN)w[1]; nu=(GEN)w[2]; if (cmpis((GEN)w[3],1)==1) return Decomppadic(p,r,fx,mf,alpha,chi,nu); } Da=lgef(nu)-3; Na=n/Da; if (mf+1<=padicprec(chi,p)) { Dchi=modii(discsr(polmodi_keep(chi,pmf)), pmf); if (gcmp0(Dchi)) Dchi=discsr(chi); } else Dchi=discsr(chi); if (gcmp0(Dchi)) alpha=gadd(alpha,gmul(p,polx[v])); else { if (vstar(p,chi)[0] > 0) alpha=gadd(alpha,gun); else { eta=setup(p,chi,polx[v],nu, &La, &Ma); if (La>1) alpha=gadd(alpha,eleval(fx,eta,alpha)); else { if (Ma==Na) return NULL; /* correspond to [fx; 1] */ w=bsrch(p,chi,ggval(Dchi,p),eta,Ma); phi=eleval(fx,(GEN)w[2],alpha); if (gcmp1((GEN)w[1])) return Decomppadic(p,r,fx,mf,phi,(GEN)w[3],(GEN)w[4]); alpha=phi; } } } } } static GEN squarefree(GEN f, GEN *ex) { GEN T,V,W,A,B; long n,i,k; T=ggcd(derivpol(f),f); V=gdeuc(f,T); n=lgef(f)-2; A=cgetg(n,t_COL); B=cgetg(n,t_COL); k=1; i=1; do { W=ggcd(T,V); T=gdeuc(T,W); if (lgef(V) != lgef(W)) { A[i]=ldeuc(V,W); B[i]=k; i++; } k++; V=W; } while (lgef(V)>3); setlg(A,i); *ex=B; return A; } GEN factorpadic4(GEN f,GEN p,long prec) { GEN w,g,poly,fx,y,p1,p2,ex,pols,exps,ppow,lead; long v=varn(f),n=lgef(f)-3,av,tetpil,mfx,i,k,j,m,r,pr; if (typ(f)!=t_POL) err(notpoler,"factorpadic"); if (gcmp0(f)) err(zeropoler,"factorpadic"); if (prec<=0) err(rootper4); if (n==0) return trivfact(); av=avma; f = padic_pol_to_int(f); if (n==1) return gerepileupto(av, padic_trivfact(f,p,prec)); f = pol_to_monic(f, &lead); pr = lead? prec + (n-1) * ggval(lead,p): prec; poly=squarefree(f,&ex); pols=cgetg(n+1,t_COL); exps=cgetg(n+1,t_COL); n = lg(poly); for (j=1,i=1; i> nsh; l=d2; while (!y[l]) l--; } while (l<=2); l++; y[1] = mymyrand() >> nsh; p1 = cgetg(BITS_IN_LONG+2,t_POL); p1[1] = evalsigne(1) | evalvarn(va); for (i=1; i>=1; } while (m); setlgef(p1,l1); y[i]=(long)to_fq(p1,a,pg); } *ptres = (GEN)y[1]; y[1] = evalsigne(1) | evalvarn(v) | evallgef(l); return y; } static void split9(GEN m, GEN *t, long d, GEN pg, GEN q, GEN munfq, GEN qq, GEN a) { long l,dv,v,av0,av,tetpil,p; GEN w,res,polmod; dv=lgef(*t)-3; if (dv==d) return; v=varn(*t); m=setloop(m); av0=avma; p = pg[2]; polmod=cgetg(3,t_POLMOD); polmod[1]=(long)dummyclone(*t); for(av=avma;;avma=av) { if (p==2) { polmod[2] = lres(stopoly92(pg,d,v,a,&res), *t); w = polmod; for (l=1; l ly) return 1; if (lx < ly) return -1; for (i=lx-1; i>1; i--) if ((fl = cmp_coeff((GEN)x[i], (GEN)y[i]))) return fl; return 0; } GEN factmod9(GEN f, GEN pp, GEN a) { long av = avma, tetpil,p,i,j,k,d,e,vf,va,nbfact,nbf,pk; GEN ex,y,f2,f3,df1,df2,g,g1,xmod,u,v,pd,q,qq,unfp,unfq,munfq,tokill, *t; GEN frobinv = gpowgs(pp, lgef(a)-4); if (typ(a)!=t_POL || typ(f)!=t_POL || gcmp0(a)) err(typeer,"factmod9"); vf=varn(f); va=varn(a); if (va<=vf) err(talker,"polynomial variable must be of higher priority than finite field\nvariable in factorff"); p=itos(pp); unfp=gmodulss(1,p); a=gmul(unfp,a); unfq=gmodulo(gmul(unfp,polun[va]), a); tokill = (GEN)unfq[1]; f = gmul(unfq,f); if (!signe(f)) err(zeropoler,"factmod9"); d = lgef(f)-3; if (!d) { avma=av; gunclone(tokill); return trivfact(); } pp = gmael(a,2,1); /* out of the stack */ t = (GEN*)cgetg(d+1,t_VEC); ex = new_chunk(d+1); xmod = cgetg(3,t_POLMOD); xmod[2] = lmul(polx[vf],unfq); munfq = gneg(unfq); qq=gpuigs(pp,lgef(a)-3); e = nbfact = 1; pk=1; df1=derivpol(f); f3=NULL; for(;;) { while (gcmp0(df1)) { pk *= p; e=pk; j=(lgef(f)-3)/p+3; setlg(f,j); setlgef(f,j); for (i=2; i>1; d++) { pd=mulii(pd,qq); v=powgi(v,qq); g=ggcd((GEN)gsub(v,xmod)[2],u); if (lgef(g) > 3) { /* Ici g est produit de pol irred ayant tous le meme degre d; */ j = nbfact+(lgef(g)-3)/d; t[nbfact]=g; q=shifti(subis(pd,1),-1); /* le premier parametre est un entier variable m qui sera * converti en un polynome w dont les coeff sont ses digits en * base p (initialement m = p --> X) pour faire pgcd de g avec * w^(p^d-1)/2 jusqu'a casser. */ split9(addis(qq,1),t+nbfact,d,pp,q,munfq,qq,a); for (; nbfact3) { t[nbfact]=u; ex[nbfact++]=e; } if (lgef(f2) == 3) break; f=f2; df1=df2; e += pk; } nbf=nbfact; tetpil=avma; y=cgetg(3,t_MAT); for (j=1; j0) g=0; } else if (ti != t_INT && ti != t_INTMOD && !is_frac_t(ti)) g=0; } l1=avma; p2=cgetg(3,t_COMPLEX); p2[1]=lmppi(DEFAULTPREC); p2[2]=ldivrs((GEN)p2[1],10); p11=cgetg(4,t_POL); p11[1]=evalsigne(1)+evallgef(4); setvarn(p11,v); p11[3]=un; p12=cgetg(5,t_POL); p12[1]=evalsigne(1)+evallgef(5); setvarn(p12,v); p12[4]=un; for (i=2; i<=deg0+2 && gcmp0((GEN)x[i]); i++) gaffsg(0,(GEN)y[i-1]); k=i-2; if (k!=deg0) { if (k) { j=deg0+3-k; pax=cgetg(j,t_POL); pax[1]=evalsigne(1)+evallgef(j); setvarn(pax,v); for (i=2; iemax) emax=e1; } e=emax-e; if (e<0) e=0; avma=l3; if (ps!=pax) xd0=deriv(ps,v); xdabs=cgetg(deg+2,t_POL); xdabs[1]=xd0[1]; for (i=2; i= -20) { p6=poleval(xd,p3); p7=gneg_i(gdiv(p4,p6)); f=1; l3=avma; if (gcmp0(p5)) exc= -32; else exc=expo(gnorm(p7))-expo(gnorm(p3)); avma=l3; for (j=1; j<=10 && f; j++) { p8=gadd(p3,p7); p9=poleval(xc,p8); p10=gnorm(p9); f=(cmprr(p10,p5)>=0)&&(exc>= -20); if (f){ gshiftz(p7,-2,p7); avma=l3; } } if (f) { avma=av1; if (DEBUGLEVEL) { fprintferr("too many iterations in rootsold(): "); fprintferr("using roots2()\n"); flusherr(); } return roots2(x,l); } else { GEN *gptr[3]; p3=p8; p4=p9; p5=p10; gptr[0]=&p3; gptr[1]=&p4; gptr[2]=&p5; gerepilemanysp(l2,l3,gptr,3); } } p1=(GEN)y[k+m*i]; setlg(p1[1],3); setlg(p1[2],3); gaffect(p3,p1); avma=l2; p14=(GEN)(p1[1]); p15=(GEN)(p1[2]); for (ln=4; ln<=l; ln=(ln<<1)-2) { setlg(p14,ln); setlg(p15,ln); if (gcmp0(p14)) { settyp(p14,t_INT); p14[1]=2; } if (gcmp0(p15)) { settyp(p15,t_INT); p15[1]=2; } p4=poleval(xc,p1); p5=poleval(xd,p1); p6=gneg_i(gdiv(p4,p5)); settyp(p14,t_REAL); settyp(p15,t_REAL); gaffect(gadd(p1,p6),p1); avma=l2; } } setlg(p14,l); setlg(p15,l); p7=gcopy(p1); p14=(GEN)(p7[1]); p15=(GEN)(p7[2]); setlg(p14,l+1); setlg(p15,l+1); if (gcmp0(p14)) { settyp(p14,t_INT); p14[1]=2; } if (gcmp0(p15)) { settyp(p15,t_INT); p15[1]=2; } for (ii=1; ii<=5; ii++) { p4=poleval(ps,p7); p5=poleval(xd0,p7); p6=gneg_i(gdiv(p4,p5)); p7=gadd(p7,p6); p14=(GEN)(p7[1]); p15=(GEN)(p7[2]); if (gcmp0(p14)) { settyp(p14,t_INT); p14[1]=2; } if (gcmp0(p15)) { settyp(p15,t_INT); p15[1]=2; } } gaffect(p7,p1); p4=poleval(ps,p7); p6=gdiv(p4,poleval(xdabs,gabs(p7,l))); if (gexpo(p6)>=expmin) { avma=av1; if (DEBUGLEVEL) { fprintferr("internal error in rootsold(): using roots2()\n"); flusherr(); } return roots2(x,l); } avma=l2; if (expo(p1[2])1) { for (j=2; j<=deg0; j++) { p1=(GEN)y[j]; if (gcmp0((GEN)p1[2])) fr=0; else fr=1; for (k=j-1; k>=1; k--) { if (gcmp0((GEN)((GEN)y[k])[2])) f=0; else f=1; if (f0) g=0; } else if (!is_const_t(ti) || ti==t_PADIC || ti==t_COMPLEX) g=0; } l1=avma; p2=cgetg(3,t_COMPLEX); p2[1]=lmppi(l); p2[2]=ldivrs((GEN)p2[1],10); p11=cgetg(4,t_POL); p11[1]=evalsigne(1)+evallgef(4); setvarn(p11,v); p11[3]=un; p12=cgetg(5,t_POL); p12[1]=evalsigne(1)+evallgef(5); setvarn(p12,v); p12[4]=un; for (i=2; (i<=deg0+2)&&(gcmp0((GEN)x[i])); i++) gaffsg(0,(GEN)y[i-1]); k=i-2; if (k!=deg0) { if (k) { j=deg0+3-k; pax=cgetg(j,t_POL); pax[1]=evalsigne(1)+evallgef(j); setvarn(pax,v); for (i=2; iemax) emax=e1; } } e=emax-e; if (e<0) e=0; avma=l3; if (ps!=pax) xd0=deriv(ps,v); xdabs=cgetg(deg+2,t_POL); xdabs[1]=xd0[1]; for (i=2; i= -20) { p6=poleval(xd,p3); p7=gneg_i(gdiv(p4,p6)); f=1; l3=avma; if (gcmp0(p5)) exc= -32; else exc=expo(gnorm(p7))-expo(gnorm(p3)); avma=l3; for (j=1; (j<=50)&&f; j++) { p8=gadd(p3,p7); p9=poleval(xc,p8); p10=gnorm(p9); f=(cmprr(p10,p5)>=0)&&(exc>= -20); if (f) { gshiftz(p7,-2,p7); avma=l3; } } if (f) err(poler9); else { GEN *gptr[3]; p3=p8; p4=p9; p5=p10; gptr[0]=&p3; gptr[1]=&p4; gptr[2]=&p5; gerepilemanysp(l2,l3,gptr,3); } } p1=(GEN)y[k+m*i]; gaffect(p3,p1); avma=l2; p14=(GEN)(p1[1]); p15=(GEN)(p1[2]); for (ln=4; ln<=l; ln=(ln<<1)-2) { if (gcmp0(p14)) { settyp(p14,t_INT); p14[1]=2; } if (gcmp0(p15)) { settyp(p15,t_INT); p15[1]=2; } p4=poleval(xc,p1); p5=poleval(xd,p1); p6=gneg_i(gdiv(p4,p5)); settyp(p14,t_REAL); settyp(p15,t_REAL); gaffect(gadd(p1,p6),p1); avma=l2; } } p7=gcopy(p1); p14=(GEN)(p7[1]); p15=(GEN)(p7[2]); setlg(p14,l+1); setlg(p15,l+1); if (gcmp0(p14)) { settyp(p14,t_INT); p14[1]=2; } if (gcmp0(p15)) { settyp(p15,t_INT); p15[1]=2; } for (ii=1; ii<=max(32,((e<=expmin)) { avma=av1; if (DEBUGLEVEL) { fprintferr("internal error in roots: using roots2\n"); flusherr(); } return roots2(x,l); } avma=l2; if ((expo(p1[2])1) { for (j=2; j<=deg0; j++) { p1=(GEN)y[j]; if (gcmp0((GEN)p1[2])) fr=0; else fr=1; for (k=j-1; k>=1; k--) { if (gcmp0((GEN)((GEN)y[k])[2])) f=0; else f=1; if (f0)? 0 : 1; } } rr=cgetg(N+1,t_COL); for (i=1; i<=N; i++) { rr[i]=lgetg(3,t_COMPLEX); p1=(GEN)rr[i]; mael(rr,i,1)=lgetr(PREC); mael(rr,i,2)=lgetr(PREC); for (j=3; j=1; j--) { x=gzero; flagrealrac=0; if (j==1) x=gneg_i(gdiv((GEN)qolbis[2],(GEN)qolbis[3])); else { x=laguer(qolbis,j,x,EPS,PREC); if (x == NULL) goto RLAB; } if (flagexactpol) { x=gprec(x,(long)((PREC-1)*pariK)); x=laguer(qol,deg,x,gmul2n(EPS,-32),PREC+1); } else x=laguer(qol,deg,x,EPS,PREC); if (x == NULL) goto RLAB; if (typ(x)==t_COMPLEX && gcmp(gabs(gimag(x),PREC),gmul2n(gmul(EPS,gabs(greal(x),PREC)),1))<=0) { x[2]=zero; flagrealrac=1; } else if (j==1 && flagrealpol) { x[2]=zero; flagrealrac=1; } else if (typ(x)!=t_COMPLEX) flagrealrac=1; for (i=1; i<=multiqol; i++) gaffect(x,(GEN)rr[nbroot+i]); nbroot+=multiqol; if (!flagrealpol || flagrealrac) { ad = (GEN*) new_chunk(j+1); for (i=0; i<=j; i++) ad[i]=(GEN)qolbis[i+2]; b=(GEN)ad[j]; for (i=j-1; i>=0; i--) { c=(GEN)ad[i]; ad[i]=b; b=gadd(gmul((GEN)rr[nbroot],b),c); } v=cgetg(j+1,t_VEC); for (i=1; i<=j; i++) v[i]=(long)ad[j-i]; qolbis=gtopoly(v,varn(qolbis)); if (flagrealpol) for (i=2; i<=j+1; i++) if (typ(qolbis[i])==t_COMPLEX) mael(qolbis,i,2)=zero; } else { ad = (GEN*) new_chunk(j-1); ad[j-2]=(GEN)qolbis[j+2]; p1=gmulsg(2,greal((GEN)rr[nbroot])); p2=gnorm((GEN)rr[nbroot]); ad[j-3]=gadd((GEN)qolbis[j+1],gmul(p1,ad[j-2])); for (i=j-2; i>=2; i--) ad[i-2] = gadd((GEN)qolbis[i+2],gsub(gmul(p1,ad[i-1]),gmul(p2,ad[i]))); v=cgetg(j,t_VEC); for (i=1; i<=j-1; i++) v[i]=(long)ad[j-1-i]; qolbis=gtopoly(v,varn(qolbis)); for (i=2; i<=j; i++) if (typ(qolbis[i])==t_COMPLEX) mael(qolbis,i,2)=zero; for (i=1; i<=multiqol; i++) gaffect(gconj((GEN)rr[nbroot]), (GEN)rr[nbroot+i]); nbroot+=multiqol; j--; } } avma=av1; } for (j=2; j<=N; j++) { x=(GEN)rr[j]; if (gcmp0((GEN)x[2])) fr=0; else fr=1; for (i=j-1; i>=1; i--) { if (gcmp0(gmael(rr,i,2))) f=0; else f=1; if (f=0; j--) { f=gadd(gmul(x,f),d); d=gadd(gmul(x,d),b); b=gadd(gmul(x,b),(GEN)pol[j+2]); erre=gadd(gnorml1(b,PREC),gmul(abx,erre)); } erre=gmul(erre,EPS); if (gcmp(gnorml1(b,PREC),erre)<=0) { gaffect(x,rac); avma = av1; return rac; } g=gdiv(d,b); g2=gsqr(g); h=gsub(g2, gmul2n(gdiv(f,b),1)); sq=gsqrt(gmulsg(N-1,gsub(gmulsg(N,h),g2)),PREC); gp=gadd(g,sq); gm=gsub(g,sq); abp=gnorm(gp); abm=gnorm(gm); if (gcmp(abp,abm)<0) gp=gcopy(gm); if (gsigne(gmax(abp,abm))==1) dx = gdivsg(N,gp); else dx = gmul(gadd(gun,abx),gexp(gmulgs(I,iter),PREC)); x1=gsub(x,dx); if (gcmp(gnorml1(gsub(x,x1),PREC),EPS)<0) { gaffect(x,rac); avma = av1; return rac; } if (iter%MT) x=gcopy(x1); else x=gsub(x,gmul(ffrac[iter/MT],dx)); } avma=av; return NULL; } #undef MR #undef MT static GEN gnorml1(GEN x,long PREC) { long av,tetpil,lx,i; GEN p1,p2,s; av=avma; switch(typ(x)) { case t_INT: case t_REAL: case t_FRAC: case t_FRACN: return gabs(x,PREC); case t_INTMOD: case t_PADIC: case t_POLMOD: case t_POL: case t_SER: case t_RFRAC: case t_RFRACN: case t_QFR: case t_QFI: return gcopy(x); case t_COMPLEX: p1=gabs((GEN)x[1],PREC); p2=gabs((GEN)x[2],PREC); tetpil=avma; return gerepile(av,tetpil,gadd(p1,p2)); case t_QUAD: p1=gabs((GEN)x[2],PREC); p2=gabs((GEN)x[3],PREC); tetpil=avma; return gerepile(av,tetpil,gadd(p1,p2)); case t_VEC: case t_COL: case t_MAT: lx=lg(x); s=gzero; for (i=1; i0){ f=gmul2n(f,-RADIX); c=gmul2n(c,-sqrdx); } if (gcmp(gdiv(gadd(c,r),f),gmul(cofgen,s))<0) { last=0; g=ginv(f); for (j=1; j<=n; j++) coeff(a,i,j)=lmul(gcoeff(a,i,j),g); for (j=1; j<=n; j++) coeff(a,j,i)=lmul(gcoeff(a,j,i),f); } } } } tetpil=avma; return gerepile(av,tetpil,gcopy(a)); } #define SIGN(a,b) ((b)>=0.0 ? fabs(a) : -fabs(a)) static GEN hqr(GEN mat) /* find all the eigenvalues of the matrix mat */ { long nn,n,m,l,k,j,its,i,mmin,flj,flk; double **a,p,q,r,s,t,u,v,w,x,y,z,anorm,*wr,*wi,eps; GEN eig; eps=0.000001; n=lg(mat)-1; a=(double**)gpmalloc(sizeof(double*)*(n+1)); for (i=1; i<=n; i++) a[i]=(double*)gpmalloc(sizeof(double)*(n+1)); for (j=1; j<=n; j++) for (i=1; i<=n; i++) a[i][j]=gtodouble((GEN)((GEN)mat[j])[i]); wr=(double*)gpmalloc(sizeof(double)*(n+1)); wi=(double*)gpmalloc(sizeof(double)*(n+1)); anorm=fabs(a[1][1]); for (i=2; i<=n; i++) for (j=(i-1); j<=n; j++) anorm+=fabs(a[i][j]); nn=n; t=0.0; if (DEBUGLEVEL>3) { fprintferr("* Finding eigenvalues\n"); flusherr(); } while (nn>=1) { its=0; do { for (l=nn; l>=2; l--) { s=fabs(a[l-1][l-1])+fabs(a[l][l]); if (s==0.0) s=anorm; if ((double)(fabs(a[l][l-1])+s)==s) break; } x=a[nn][nn]; if (l==nn){ wr[nn]=x+t; wi[nn--]=0.0; } else { y=a[nn-1][nn-1]; w=a[nn][nn-1]*a[nn-1][nn]; if (l==(nn-1)) { p=0.5*(y-x); q=p*p+w; z=sqrt(fabs(q)); x+=t; if ((q>=0.0)||(fabs(q)<=eps)) { z=p+SIGN(z,p); wr[nn-1]=wr[nn]=x+z; if (fabs(z)>eps) wr[nn]=x-w/z; wi[nn-1]=wi[nn]=0.0; } else{ wr[nn-1]=wr[nn]=x+p; wi[nn-1]=-(wi[nn]=z); } nn-=2; } else { if (its==30) err(talker,"too many iterations in hqr"); if ((its==10)||(its==20)) { t+=x; for (i=1; i<=nn; i++) a[i][i]-=x; s=fabs(a[nn][nn-1])+fabs(a[nn-1][nn-2]); y=x=0.75*s; w=-0.4375*s*s; } ++its; for (m=(nn-2); m>=l; m--) { z=a[m][m]; r=x-z; s=y-z; p=(r*s-w)/a[m+1][m]+a[m][m+1]; q=a[m+1][m+1]-z-r-s; r=a[m+2][m+1]; s=fabs(p)+fabs(q)+fabs(r); p/=s; q/=s; r/=s; if (m==l) break; u=fabs(a[m][m-1])*(fabs(q)+fabs(r)); v=fabs(p)*(fabs(a[m-1][m-1])+fabs(z)+fabs(a[m+1][m+1])); if ((double)(u+v)==v) break; } for (i=m+2; i<=nn; i++){ a[i][i-2]=0.0; if (i!=(m+2)) a[i][i-3]=0.0; } for (k=m; k<=nn-1; k++) { if (k!=m) { p=a[k][k-1]; q=a[k+1][k-1]; r=0.0; if (k!=(nn-1)) r=a[k+2][k-1]; if ((x=fabs(p)+fabs(q)+fabs(r))!=0.0){ p/=x; q/=x; r/=x; } } if ((s=SIGN(sqrt(p*p+q*q+r*r),p))!=0.0) { if (k==m){ if (l!=m) a[k][k-1]=-a[k][k-1]; }else a[k][k-1]=-s*x; p+=s; x=p/s; y=q/s; z=r/s; q/=p; r/=p; for (j=k; j<=nn; j++) { p=a[k][j]+q*a[k+1][j]; if (k!=(nn-1)){ p+=r*a[k+2][j]; a[k+2][j]-=p*z; } a[k+1][j]-=p*y; a[k][j]-=p*x; } mmin=(nn=1; k--) { if (wi[k]) flk=1; else flk=0; if (flk0)) break; wr[k+1]=wr[k]; wi[k+1]=wi[k]; } wr[k+1]=x; wi[k+1]=y; } if (DEBUGLEVEL>3) { fprintferr("* End of the computation of eigenvalues\n"); flusherr(); } for (i=1; i<=n; i++) free(a[i]); free(a); eig=cgetg(n+1,t_COL); for (i=1; i<=n; i++) { if (wi[i]) { eig[i]=lgetg(3,t_COMPLEX); ((GEN)eig[i])[1]=(long)dbltor(wr[i]); ((GEN)eig[i])[2]=(long)dbltor(wi[i]); } else eig[i]=(long)dbltor(wr[i]); } free(wr); free(wi); return eig; } static GEN zrhqr(GEN a,long PREC) /* ONLY FOR POLYNOMIAL WITH REAL COEFFICIENTS : give the roots of * the polynomial a (first, the real roots, then the * non real roots) in increasing order of their real * parts MULTIPLE ROOTS ARE FORBIDDEN. */ { long av,tetpil,n,i,j,k,ct,prec; GEN aa,b,p1,rt,rr,hess,x,dx,y,hessbis,eps,newval; GEN oldval = NULL; /* for lint */ av=avma; n=lgef(a)-3; prec=PREC; hess=cgetg(n+1,t_MAT); for (k=1; k<=n; k++) hess[k]=lgetg(n+1,t_COL); for (k=1; k<=n; k++) { p1=(GEN)hess[k]; p1[1]=lneg(gdiv((GEN)a[n-k+2],(GEN)a[n+2])); for (j=2; j<=n; j++){ if (j==(k+1)) p1[j]=un; else p1[j]=zero; } } rr=cgetg(n+1,t_COL); for (i=1; i<=n; i++) { rr[i]=lgetg(3,t_COMPLEX); mael(rr,i,1)=lgetr(PREC); mael(rr,i,2)=lgetr(PREC); } hessbis=balanc(hess); rt=hqr(hessbis); eps=cgetr(prec); p1=gpuigs(gdeux,-bit_accuracy(prec)); gaffect(p1,eps); prec=2*PREC; /* polishing the roots */ aa=cgetg(n+3,t_POL); aa[1]=a[1]; for (i=2; i<=n+2; i++){ aa[i]=lgetr(prec); gaffect((GEN)a[i],(GEN)aa[i]); } b=deriv(aa,varn(aa)); for (i=1; i<=n; i++) { ct=0; if (typ(rt[i])==t_REAL) { x=cgetr(prec); affrr((GEN)rt[i],x); } else { x=cgetg(3,t_COMPLEX); x[1]=lgetr(prec); affrr(gmael(rt,i,1),(GEN)x[1]); x[2]=lgetr(prec); affrr(gmael(rt,i,2),(GEN)x[2]); } LAB1: dx=poleval(b,x); if (gcmp(gabs(dx,prec),eps) <= 0) err(talker,"the polynomial has probably multiple roots in zrhqr"); y=gsub(x,gdiv(poleval(aa,x),dx)); newval=gabs(poleval(aa,y),prec); if (gcmp(newval,eps) > 0 && (!ct || gcmp(newval,oldval) < 0)) { ct++; oldval=newval; x=y; goto LAB1; } gaffect(y,(GEN)rr[i]); } if (DEBUGLEVEL>3){ fprintferr("polished roots = %Z",rr); flusherr(); } tetpil=avma; return gerepile(av,tetpil,gcopy(rr)); }