/* $Id: gen3.c,v 1.39 2001/10/01 12:11:31 karim Exp $ Copyright (C) 2000 The PARI group. This file is part of the PARI/GP package. PARI/GP is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation. It is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY WHATSOEVER. Check the License for details. You should have received a copy of it, along with the package; see the file 'COPYING'. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. */ /********************************************************************/ /** **/ /** GENERIC OPERATIONS **/ /** (third part) **/ /** **/ /********************************************************************/ #include "pari.h" /********************************************************************/ /** **/ /** PRINCIPAL VARIABLE NUMBER **/ /** **/ /********************************************************************/ int gvar(GEN x) { long tx=typ(x),i,v,w; if (is_polser_t(tx)) return varn(x); if (tx==t_POLMOD) return varn(x[1]); if (is_const_t(tx) || is_qf_t(tx) || is_noncalc_t(tx)) return BIGINT; v=BIGINT; for (i=1; i < lg(x); i++) { w=gvar((GEN)x[i]); if (w>TWOPOTBITS_IN_LONG); l=lg(x); if (l>k) k=l; return k; } if (tx==t_COMPLEX) { k=precision((GEN)x[1]); l=precision((GEN)x[2]); if (l && (!k || l=0, par rapport a la variable * principale si v<0. On suppose x fraction rationnelle ou polynome. * Convention deg(0)=-1. */ long poldegree(GEN x, long v) { long tx=typ(x), av, w, d; if (is_scalar_t(tx)) return gcmp0(x)? -1: 0; switch(tx) { case t_POL: w = varn(x); if (v < 0 || v == w) return degpol(x); if (v < w) return signe(x)? 0: -1; av = avma; x = gsubst(gsubst(x,w,polx[MAXVARN]),v,polx[0]); if (gvar(x)) { d = gcmp0(x)? -1: 0; } else d = degpol(x); avma = av; return d; case t_RFRAC: case t_RFRACN: if (gcmp0((GEN)x[1])) return -1; return poldegree((GEN)x[1],v) - poldegree((GEN)x[2],v); } err(typeer,"degree"); return 0; /* not reached */ } long degree(GEN x) { return poldegree(x,-1); } /* si v<0, par rapport a la variable principale, sinon par rapport a v. * On suppose que x est un polynome ou une serie. */ GEN pollead(GEN x, long v) { long l,tx = typ(x),av,tetpil,w; GEN xinit; if (is_scalar_t(tx)) return gcopy(x); w = varn(x); switch(tx) { case t_POL: if (v < 0 || v == w) { l=lgef(x); return (l==2)? gzero: gcopy((GEN)x[l-1]); } if (v < w) return gcopy(x); av = avma; xinit = x; x = gsubst(gsubst(x,w,polx[MAXVARN]),v,polx[0]); if (gvar(x)) { avma = av; return gcopy(xinit); } l = lgef(x); if (l == 2) { avma = av; return gzero; } tetpil = avma; x = gsubst((GEN)x[l-1],MAXVARN,polx[w]); return gerepile(av,tetpil,x); case t_SER: if (v < 0 || v == w) return (!signe(x))? gzero: gcopy((GEN)x[2]); if (v < w) return gcopy(x); av = avma; xinit = x; x = gsubst(gsubst(x,w,polx[MAXVARN]),v,polx[0]); if (gvar(x)) { avma = av; return gcopy(xinit);} if (!signe(x)) { avma = av; return gzero;} tetpil = avma; x = gsubst((GEN)x[2],MAXVARN,polx[w]); return gerepile(av,tetpil,x); } err(typeer,"pollead"); return NULL; /* not reached */ } /* returns 1 if there's a real component in the structure, 0 otherwise */ int isinexactreal(GEN x) { long tx=typ(x),lx,i; if (is_scalar_t(tx)) { switch(tx) { case t_REAL: return 1; case t_COMPLEX: case t_QUAD: return (typ(x[1])==t_REAL || typ(x[2])==t_REAL); } return 0; } switch(tx) { case t_QFR: case t_QFI: return 0; case t_RFRAC: case t_RFRACN: return isinexactreal((GEN)x[1]) || isinexactreal((GEN)x[2]); } if (is_noncalc_t(tx)) return 0; lx = (tx==t_POL)? lgef(x): lg(x); for (i=lontyp[tx]; i1; i--) if (!isexactzero((GEN)g[i])) return 0; return 1; case t_VEC: case t_COL: case t_MAT: for (i=lg(g)-1; i; i--) if (!isexactzero((GEN)g[i])) return 0; return 1; } return 0; } int iscomplex(GEN x) { switch(typ(x)) { case t_INT: case t_REAL: case t_FRAC: case t_FRACN: return 0; case t_COMPLEX: return !gcmp0((GEN)x[2]); case t_QUAD: return signe(mael(x,1,2)) > 0; } err(typeer,"iscomplex"); return 0; /* not reached */ } int ismonome(GEN x) { long i; if (typ(x)!=t_POL || !signe(x)) return 0; for (i=lgef(x)-2; i>1; i--) if (!isexactzero((GEN)x[i])) return 0; return 1; } /********************************************************************/ /** **/ /** MULTIPLICATION SIMPLE **/ /** **/ /********************************************************************/ #define fix_frac(z) if (signe(z[2])<0)\ {\ setsigne(z[1],-signe(z[1]));\ setsigne(z[2],1);\ } /* assume z[1] was created last */ #define fix_frac_if_int(z) if (is_pm1(z[2]))\ z = gerepileuptoint((long)(z+3), (GEN)z[1]); GEN gmulsg(long s, GEN y) { long ty=typ(y),ly=lg(y),i,av,tetpil; GEN z,p1,p2; switch(ty) { case t_INT: return mulsi(s,y); case t_REAL: return mulsr(s,y); case t_INTMOD: z=cgetg(3,t_INTMOD); p2=(GEN)y[1]; (void)new_chunk(lgefint(p2)<<2); /* HACK */ p1=mulsi(s,(GEN)y[2]); avma=(long)z; z[2]=lmodii(p1,p2); icopyifstack(p2,z[1]); return z; case t_FRAC: if (!s) return gzero; z = cgetg(3,t_FRAC); i = cgcd(s, smodis((GEN)y[2], s)); if (i == 1) { z[2] = licopy((GEN)y[2]); z[1] = lmulis((GEN)y[1], s); } else { z[2] = ldivis((GEN)y[2], i); z[1] = lmulis((GEN)y[1], s/i); fix_frac_if_int(z); } return z; case t_FRACN: z=cgetg(3,ty); z[1]=lmulsi(s,(GEN)y[1]); z[2]=licopy((GEN)y[2]); return z; case t_COMPLEX: z=cgetg(ly,ty); z[1]=lmulsg(s,(GEN)y[1]); z[2]=lmulsg(s,(GEN)y[2]); return z; case t_PADIC: if (!s) return gzero; av=avma; p1=cgetp(y); gaffsg(s,p1); tetpil=avma; return gerepile(av,tetpil,gmul(p1,y)); case t_QUAD: z=cgetg(ly,ty); copyifstack(y[1],z[1]); z[2]=lmulsg(s,(GEN)y[2]); z[3]=lmulsg(s,(GEN)y[3]); return z; case t_POLMOD: z=cgetg(ly,ty); z[2]=lmulsg(s,(GEN)y[2]); copyifstack(y[1],z[1]); return z; case t_POL: if (!s || !signe(y)) return zeropol(varn(y)); ly=lgef(y); z=cgetg(ly,t_POL); z[1]=y[1]; for (i=2; i varn(y)) return (lgef(y) < 4)? gzero: gcopy(x); if (varn(x[1])!=varn(y)) return gzero; z=cgetg(3,t_POLMOD); z[1]=lgcd((GEN)x[1],y); z[2]=lres((GEN)x[2],(GEN)z[1]); return z; } switch(tx) { case t_POL: return gres(x,y); case t_RFRAC: case t_RFRACN: av=avma; if (tx==t_RFRACN) x=gred_rfrac(x); p1=gmul((GEN)x[1],ginvmod((GEN)x[2],y)); tetpil=avma; return gerepile(av,tetpil,gres(p1,y)); case t_SER: if (ismonome(y) && varn(x) == varn(y)) { long d = degpol(y); if (lg(x)-2 + valp(x) < d) err(operi,"%",tx,ty); av = avma; return gerepileupto(av, gmod(gtrunc(x), y)); } default: err(operf,"%",tx,ty); } } err(operf,"%",tx,ty); return NULL; /* not reached */ } GEN gmodulsg(long x, GEN y) { GEN z; switch(typ(y)) { case t_INT: z=cgetg(3,t_INTMOD); z[1]=labsi(y); z[2]=lmodsi(x,y); return z; case t_POL: z=cgetg(3,t_POLMOD); z[1]=lcopy(y); z[2]=lstoi(x); return z; } err(operf,"%",t_INT,typ(y)); return NULL; /* not reached */ } GEN gmodulss(long x, long y) { GEN z=cgetg(3,t_INTMOD); y=labs(y); z[1]=lstoi(y); z[2]=lstoi(x % y); return z; } static GEN specialmod(GEN x, GEN y) { GEN z = gmod(x,y); if (gvar(z) < varn(y)) z = gzero; return z; } GEN gmodulo(GEN x,GEN y) { long tx=typ(x),l,i; GEN z; if (is_matvec_t(tx)) { l=lg(x); z=cgetg(l,tx); for (i=1; i 3)? lcopy(x): lmod(x,y); return z; } if (tx!=t_POL && !is_rfrac_t(tx) && tx!=t_SER) break; z[2]=(long)specialmod(x,y); return z; } err(operf,"%",tx,typ(y)); return NULL; /* not reached */ } GEN Mod0(GEN x,GEN y,long flag) { switch(flag) { case 0: return gmodulcp(x,y); case 1: return gmodulo(x,y); default: err(flagerr,"Mod"); } return NULL; /* not reached */ } /*******************************************************************/ /* */ /* DIVISION ENTIERE GENERALE */ /* DIVISION ENTIERE AVEC RESTE GENERALE */ /* */ /*******************************************************************/ GEN gdivent(GEN x, GEN y) { long tx=typ(x),ty=typ(y); if (tx == t_INT) { if (ty == t_INT) return truedvmdii(x,y,NULL); if (ty!=t_POL) err(typeer,"gdivent"); return gzero; } if (ty != t_POL) err(typeer,"gdivent"); if (tx == t_POL) return gdeuc(x,y); if (! is_scalar_t(tx)) err(typeer,"gdivent"); return gzero; } GEN gdiventres(GEN x, GEN y) { long tx=typ(x),ty=typ(y); GEN z = cgetg(3,t_COL); if (tx==t_INT) { if (ty==t_INT) { z[1]=(long)truedvmdii(x,y,(GEN*)(z+2)); return z; } if (ty!=t_POL) err(typeer,"gdiventres"); z[1]=zero; z[2]=licopy(x); return z; } if (ty != t_POL) err(typeer,"gdiventres"); if (tx == t_POL) { z[1]=ldivres(x,y,(GEN *)(z+2)); return z; } if (! is_scalar_t(tx)) err(typeer,"gdiventres"); z[1]=zero; z[2]=lcopy(x); return z; } GEN gdivmod(GEN x, GEN y, GEN *pr) { long ty,tx=typ(x); if (tx==t_INT) { ty=typ(y); if (ty==t_INT) return dvmdii(x,y,pr); if (ty==t_POL) { *pr=gcopy(x); return gzero; } err(typeer,"gdivmod"); } if (tx!=t_POL) err(typeer,"gdivmod"); return poldivres(x,y,pr); } /* When x and y are integers, compute the quotient x/y, rounded to the * nearest integer. If there is a tie, the quotient is rounded towards zero. * If x and y are not both integers, same as gdivent. */ GEN gdivround(GEN x, GEN y) { long tx=typ(x),ty=typ(y); if (tx==t_INT) { if (ty==t_INT) { long av = avma, av1,fl; GEN r, q=dvmdii(x,y,&r); /* q = x/y rounded towards 0, sqn(r)=sgn(x) */ if (r==gzero) return q; av1 = avma; fl = absi_cmp(shifti(r,1),y); avma = av1; cgiv(r); if (fl >= 0) /* If 2*|r| >= |y| */ { long sz = signe(x)*signe(y); if (fl || sz > 0) { av1=avma; q = gerepile(av,av1,addis(q,sz)); } } return q; } if (ty!=t_POL) err(typeer,"gdivround"); return gzero; } if (ty != t_POL) err(typeer,"gdivround"); if (tx == t_POL) return gdeuc(x,y); if (! is_scalar_t(tx)) err(typeer,"gdivround"); return gzero; } /*******************************************************************/ /* */ /* SHIFT D'UN GEN */ /* */ /*******************************************************************/ /* Shift tronque si n<0 (multiplication tronquee par 2^n) */ GEN gshift(GEN x, long n) { long i,l,lx, tx = typ(x); GEN y; switch(tx) { case t_INT: return shifti(x,n); case t_REAL: return shiftr(x,n); case t_VEC: case t_COL: case t_MAT: lx=lg(x); y=cgetg(lx,tx); l=lontyp[tx]; for (i=1; i=0) return shifti(x,n); if (!signe(x)) return gzero; l = vali(x); n = -n; if (n<=l) return shifti(x,-n); y=cgetg(3,t_FRAC); y[1]=lshifti(x,-l); y[2]=lshifti(gun,n-l); return y; case t_REAL: return shiftr(x,n); case t_INTMOD: if (n > 0) { y=cgetg(3,t_INTMOD); p2=(GEN)x[1]; av=avma; new_chunk(2 * lgefint(p2) + n); /* HACK */ p1 = shifti((GEN)x[2],n); avma=av; y[2]=lmodii(p1,p2); icopyifstack(p2,y[1]); return y; } l=avma; y=gmul2n(gun,n); tetpil=avma; return gerepile(l,tetpil,gmul(y,x)); case t_FRAC: case t_FRACN: l = vali((GEN)x[1]); k = vali((GEN)x[2]); if (n+l-k>=0) { if (expi((GEN)x[2]) == k) /* x[2] power of 2 */ return shifti((GEN)x[1],n-k); l = n-k; k = -k; } else { k = -l-n; l = -l; } y=cgetg(3,t_FRAC); y[1]=lshifti((GEN)x[1],l); y[2]=lshifti((GEN)x[2],k); return y; case t_QUAD: y=cgetg(4,t_QUAD); copyifstack(x[1],y[1]); for (i=2; i<4; i++) y[i]=lmul2n((GEN)x[i],n); return y; case t_POLMOD: y=cgetg(3,t_POLMOD); copyifstack(x[1],y[1]); y[2]=lmul2n((GEN)x[2],n); return y; case t_POL: case t_COMPLEX: case t_SER: case t_VEC: case t_COL: case t_MAT: lx = (tx==t_POL)? lgef(x): lg(x); y=cgetg(lx,tx); l=lontyp[tx]; for (i=1; i=0) { y[1]=lmul2n((GEN)x[1],n); y[2]=lcopy((GEN)x[2]); } else { y[2]=lmul2n((GEN)x[2],-n); y[1]=lcopy((GEN)x[1]); } return y; case t_PADIC: l=avma; y=gmul2n(gun,n); tetpil=avma; return gerepile(l,tetpil,gmul(y,x)); } err(typeer,"gmul2n"); return NULL; /* not reached */ } /*******************************************************************/ /* */ /* INVERSE D'UN GEN */ /* */ /*******************************************************************/ extern GEN fix_rfrac_if_pol(GEN x, GEN y); GEN ginv(GEN x) { long tx=typ(x),av,tetpil,s; GEN z,y,p1,p2; switch(tx) { case t_INT: if (is_pm1(x)) return icopy(x); if (!signe(x)) err(gdiver2); z=cgetg(3,t_FRAC); z[1] = (signe(x)<0)? lnegi(gun): un; z[2]=labsi(x); return z; case t_REAL: return divsr(1,x); case t_INTMOD: z=cgetg(3,t_INTMOD); icopyifstack(x[1],z[1]); z[2]=lmpinvmod((GEN)x[2],(GEN)x[1]); return z; case t_FRAC: case t_FRACN: s = signe(x[1]); if (!s) err(gdiver2); if (is_pm1(x[1])) return s>0? icopy((GEN)x[2]): negi((GEN)x[2]); z = cgetg(3,tx); z[1] = licopy((GEN)x[2]); z[2] = licopy((GEN)x[1]); if (s < 0) { setsigne(z[1],-signe(z[1])); setsigne(z[2],1); } return z; case t_COMPLEX: case t_QUAD: av=avma; p1=gnorm(x); p2=gconj(x); tetpil=avma; return gerepile(av,tetpil,gdiv(p2,p1)); case t_PADIC: z = cgetg(5,t_PADIC); if (!signe(x[4])) err(gdiver2); z[1] = evalprecp(precp(x)) | evalvalp(-valp(x)); icopyifstack(x[2], z[2]); z[3] = licopy((GEN)x[3]); z[4] = lmpinvmod((GEN)x[4],(GEN)z[3]); return z; case t_POLMOD: z=cgetg(3,t_POLMOD); copyifstack(x[1],z[1]); z[2]=linvmod((GEN)x[2],(GEN)x[1]); return z; case t_POL: case t_SER: return gdiv(gun,x); case t_RFRAC: case t_RFRACN: if (gcmp0((GEN) x[1])) err(gdiver2); p1 = fix_rfrac_if_pol((GEN)x[2],(GEN)x[1]); if (p1) return p1; z=cgetg(3,tx); z[1]=lcopy((GEN)x[2]); z[2]=lcopy((GEN)x[1]); return z; case t_QFR: { long k,l; l=signe(x[2]); setsigne(x[2],-l); k=signe(x[4]); setsigne(x[4],-k); z=redreal(x); setsigne(x[2],l); setsigne(x[4],k); return z; } case t_QFI: y=gcopy(x); if (!egalii((GEN)x[1],(GEN)x[2]) && !egalii((GEN)x[1],(GEN)x[3])) setsigne(y[2],-signe(y[2])); return y; case t_MAT: return (lg(x)==1)? cgetg(1,t_MAT): invmat(x); } err(typeer,"inverse"); return NULL; /* not reached */ } /*******************************************************************/ /* */ /* SUBSTITUTION DANS UN POLYNOME OU UNE SERIE */ /* */ /*******************************************************************/ /* Convert t_SER --> t_POL / t_RFRAC */ static GEN gconvsp(GEN x, int flpile) { long v = varn(x), av,tetpil,i; GEN p1,y; if (gcmp0(x)) return zeropol(v); av=avma; y=dummycopy(x); settyp(y,t_POL); i=lg(x)-1; while (i>1 && gcmp0((GEN)y[i])) i--; setlgef(y,i+1); p1=gpuigs(polx[v],valp(x)); tetpil=avma; p1=gmul(p1,y); return flpile? gerepile(av,tetpil,p1): p1; } GEN gsubst0(GEN x, GEN T, GEN y) { ulong av; long d, v; if (typ(T) != t_POL || !ismonome(T) || !gcmp1(leading_term(T))) err(talker,"variable number expected in subst"); d = degpol(T); v = varn(T); if (d == 1) return gsubst(x, v, y); av = avma; return gerepilecopy(av, gsubst(gdeflate(x, v, d), v, y)); } GEN gsubst(GEN x, long v, GEN y) { long tx = typ(x), ty = typ(y), lx = lg(x), ly = lg(y); long l,vx,vy,e,ex,ey,tetpil,av,i,j,k,jb,limite; GEN t,p1,p2,z; if (ty==t_MAT) { if (ly==1) return cgetg(1,t_MAT); if (ly != lg(y[1])) err(talker,"forbidden substitution by a non square matrix"); } else if (is_graphicvec_t(ty)) err(talker,"forbidden substitution by a vector"); if (is_scalar_t(tx)) { if (tx!=t_POLMOD || v <= varn(x[1])) { if (ty==t_MAT) return gscalmat(x,ly-1); return gcopy(x); } av=avma; p1=gsubst((GEN)x[1],v,y); vx=varn(p1); p2=gsubst((GEN)x[2],v,y); vy=gvar(p2); if (typ(p1)!=t_POL) err(talker,"forbidden substitution in a scalar type"); if (vy>=vx) { tetpil=avma; return gerepile(av,tetpil,gmodulcp(p2,p1)); } lx=lgef(p2); tetpil=avma; z=cgetg(lx,t_POL); z[1]=p2[1]; for (i=2; i2; i--) z=gadd(gmul(z,p1),gsubst((GEN)x[i-1],v,y)); return gerepileupto(av,z); } if (ty!=t_MAT) return (vx>v)? gcopy(x): poleval(x,y); if (vx>v) return gscalmat(x,ly-1); if (l==3) return gscalmat((GEN)x[2],ly-1); av=avma; z=(GEN)x[l-1]; for (i=l-1; i>2; i--) z=gaddmat((GEN)x[i-1],gmul(z,y)); return gerepileupto(av,z); case t_SER: vx=varn(x); if (vx > v) return (ty==t_MAT)? gscalmat(x,ly-1): gcopy(x); ex = valp(x); if (vx < v) { if (!signe(x)) return gcopy(x); /* a ameliorer */ av=avma; setvalp(x,0); p1=gconvsp(x,0); setvalp(x,ex); p2=gsubst(p1,v,y); tetpil=avma; z=tayl(p2,vx,lx-2); if (ex) { p1=gpuigs(polx[vx],ex); tetpil=avma; z=gmul(z,p1); } return gerepile(av,tetpil,z); } switch(ty) /* here vx == v */ { case t_SER: ey=valp(y); vy=varn(y); if (ey<1) return zeroser(vy,ey*(ex+lx-2)); l=(lx-2)*ey+2; if (ex) { if (l>ly) l=ly; } else { if (gcmp0(y)) l=ey+2; else { if (l>ly) l=ey+ly; } } if (vy!=vx) { av=avma; z = zeroser(vy,0); for (i=lx-1; i>=2; i--) z = gadd((GEN)x[i], gmul(y,z)); if (ex) z = gmul(z, gpuigs(y,ex)); return gerepileupto(av,z); } av=avma; limite=stack_lim(av,1); t=cgetg(ly,t_SER); z = scalarser((GEN)x[2],varn(y),l-2); for (i=2; i1; j--) { p1=gzero; for (k=2; k1) err(warnmem,"gsubst"); gptr[0]=&z; gptr[1]=&t; gerepilemany(av,gptr,2); } } if (!ex) return gerepilecopy(av,z); if (l2; k--) p1=gsub(p1,gmul((GEN)u[k],(GEN)x[j-k+2])); u[j]=lsub(p1,(GEN)x[j]); } p1=gmulsg(i,(GEN)x[i+1]); for (k=2; k1) err(warnmem,"recip"); for(k=i+1; kv) return gzero; if (vxv) return gzero; if (vxvarn(x[1])) { y=cgetg(3,t_POLMOD); copyifstack(x[1],y[1]); y[2]=linteg((GEN)x[2],v); return y; } if (gcmp0(x)) return gzero; y=cgetg(4,t_POL); y[1] = evalsigne(1) | evallgef(4) | evalvarn(v); y[2]=zero; y[3]=lcopy(x); return y; } switch(tx) { case t_POL: vx=varn(x); lx=lgef(x); if (lx==2) return zeropol(min(v,vx)); if (vx>v) { y=cgetg(4,t_POL); y[1] = signe(x)? evallgef(4) | evalvarn(v) | evalsigne(1) : evallgef(4) | evalvarn(v); y[2]=zero; y[3]=lcopy(x); return y; } if (vxv) { y=cgetg(4,t_POL); y[1] = evallgef(4) | evalvarn(v) | evalsigne(1); y[2]=zero; y[3]=lcopy(x); return y; } if (vxv) { y=cgetg(4,t_POL); y[1] = signe(x[1])? evallgef(4) | evalvarn(v) | evalsigne(1) : evallgef(4) | evalvarn(v); y[2]=zero; y[3]=lcopy(x); return y; } if (vx0) { cgiv(p1); tetpil=avma; return gerepile(av,tetpil,addsi(1,y)); } return y; case t_RFRAC: case t_RFRACN: return gdeuc((GEN)x[1],(GEN)x[2]); case t_VEC: case t_COL: case t_MAT: lx=lg(x); y=cgetg(lx,tx); for (i=1; i0? gun: negi(gun); av=avma; p1 = realun(3 + (ex>>TWOPOTBITS_IN_LONG)); setexpo(p1,-1); /* p1 = 0.5 */ p1 = addrr(x,p1); tetpil = avma; return gerepile(av,tetpil,mpent(p1)); } case t_FRAC: case t_FRACN: av=avma; p1 = addii(shifti((GEN)x[2], -1), (GEN)x[1]); return gerepileuptoint(av, truedvmdii(p1, (GEN)x[2], NULL)); case t_POLMOD: y=cgetg(3,t_POLMOD); copyifstack(x[1],y[1]); y[2]=lround((GEN)x[2]); return y; case t_COMPLEX: case t_POL: case t_SER: case t_RFRAC: case t_RFRACN: case t_VEC: case t_COL: case t_MAT: lx = (tx==t_POL)? lgef(x): lg(x); y=cgetg(lx,tx); for (i=1; i=0) { *e=expo(x); avma=av; return gzero; } *e=expo(addsr(1,x)); avma=av; return negi(gun); } lx=lg(x); e1 = ex - bit_accuracy(lx) + 1; settyp(p1,t_INT); setlgefint(p1,lx); y=shifti(p1,e1); if (signe(x)<0) y=addsi(-1,y); y = gerepileupto(av,y); if (e1<=0) { av=avma; e1=expo(subri(x,y)); avma=av; } *e=e1; return y; case t_POLMOD: y=cgetg(3,t_POLMOD); copyifstack(x[1],y[1]); y[2]=lrndtoi((GEN)x[2],e); return y; case t_COMPLEX: case t_POL: case t_SER: case t_RFRAC: case t_RFRACN: case t_VEC: case t_COL: case t_MAT: lx=(tx==t_POL)? lgef(x): lg(x); y=cgetg(lx,tx); for (i=1; i*e) *e=e1; } return y; } err(typeer,"grndtoi"); return NULL; /* not reached */ } /* e = number of error bits on integral part */ GEN gcvtoi(GEN x, long *e) { long tx=typ(x), lx,i,ex,av,e1; GEN y; *e = -HIGHEXPOBIT; if (tx == t_REAL) { long x0,x1; ex=expo(x); if (ex<0) { *e=ex; return gzero; } lx=lg(x); e1 = ex - bit_accuracy(lx) + 1; x0=x[0]; x1=x[1]; settyp(x,t_INT); setlgefint(x,lx); y=shifti(x,e1); x[0]=x0; x[1]=x1; if (e1<=0) { av=avma; e1=expo(subri(x,y)); avma=av; } *e=e1; return y; } if (is_matvec_t(tx)) { lx=lg(x); y=cgetg(lx,tx); for (i=1; i*e) *e=e1; } return y; } return gtrunc(x); } /* smallest integer greater than any incarnations of the real x * [avoid mpent() and "precision loss in truncation"] */ GEN ceil_safe(GEN x) { ulong av = avma; long e; GEN y = gcvtoi(x,&e); if (e < 0) e = 0; y = addii(y, shifti(gun,e)); return gerepileuptoint(av, y); } GEN gtrunc(GEN x) { long tx=typ(x),av,tetpil,i,v; GEN y; switch(tx) { case t_INT: case t_POL: return gcopy(x); case t_REAL: return mptrunc(x); case t_FRAC: case t_FRACN: return divii((GEN)x[1],(GEN)x[2]); case t_PADIC: if (!signe(x[4])) return gzero; v = valp(x); if (!v) return gcopy((GEN)x[4]); if (v>0) { /* here p^v is an integer */ av=avma; y=gpuigs((GEN)x[2],v); tetpil=avma; return gerepile(av,tetpil, mulii(y,(GEN)x[4])); } y=cgetg(3,t_FRAC); y[1]=licopy((GEN)x[4]); y[2]=lpuigs((GEN)x[2],-v); return y; case t_RFRAC: case t_RFRACN: return gdeuc((GEN)x[1],(GEN)x[2]); case t_SER: return gconvsp(x,1); case t_VEC: case t_COL: case t_MAT: { long lx=lg(x); y=cgetg(lx,tx); for (i=1; i INT, POL & SER */ /* */ /*******************************************************************/ /* return a_n B^n + ... + a_0, where B = 2^BIL. Assume n even if BIL=64 and * the a_i are 32bits integers, one of which is non-zero */ GEN coefs_to_int(long n, ...) { va_list ap; GEN x, y; long i; va_start(ap,n); #ifdef LONG_IS_64BIT n >>= 1; #endif x = cgetg(n+2, t_INT); y = x + 2; x[1] = evallgefint(n+2) | evalsigne(1); for (i=0; i = 0; i--) y[i] = (long) va_arg(ap, GEN); return normalizepol(x); } GEN zeropol(long v) { GEN x = cgetg(2,t_POL); x[1] = evallgef(2) | evalvarn(v); return x; } GEN scalarpol(GEN x, long v) { GEN y=cgetg(3,t_POL); y[1] = gcmp0(x)? evallgef(3) | evalvarn(v) : evallgef(3) | evalvarn(v) | evalsigne(1); y[2]=lcopy(x); return y; } /* deg1pol(a,b,x)=a*x+b*/ GEN deg1pol(GEN x1, GEN x0,long v) { GEN x = cgetg(4,t_POL); x[1] = evalsigne(1) | evalvarn(v) | evallgef(4); x[2] = lcopy(x0); x[3] = lcopy(x1); return normalizepol_i(x,4); } /* same, no copy */ GEN deg1pol_i(GEN x1, GEN x0,long v) { GEN x = cgetg(4,t_POL); x[1] = evalsigne(1) | evalvarn(v) | evallgef(4); x[2] = (long)x0; x[3] = (long)x1; return normalizepol_i(x,4); } static GEN gtopoly0(GEN x, long v, int reverse) { long tx=typ(x),lx,i,j; GEN y; if (v<0) v = 0; if (isexactzero(x)) return zeropol(v); if (is_scalar_t(tx)) return scalarpol(x,v); switch(tx) { case t_POL: y=gcopy(x); break; case t_SER: y=gconvsp(x,1); if (typ(y) != t_POL) err(talker,"t_SER with negative valuation in gtopoly"); break; case t_RFRAC: case t_RFRACN: y=gdeuc((GEN)x[1],(GEN)x[2]); break; case t_QFR: case t_QFI: case t_VEC: case t_COL: case t_MAT: lx=lg(x); if (reverse) { while (lx-- && isexactzero((GEN)x[lx])); i=lx+2; y=cgetg(i,t_POL); y[1]=evallgef(i); if (!gcmp0(x)) y[1] |= evalsigne(1); for (j=2; jl) l=precdl; y=cgetg(l+2,t_SER); y[1] = evalsigne(1) | evalvalp(i-2) | evalvarn(v); for (j=2; j<=lx-i+1; j++) y[j]=lcopy((GEN)x[j+i-2]); for ( ; j<=l+1; j++) y[j]=zero; break; case t_RFRAC: case t_RFRACN: av=avma; p1=gtoser((GEN)x[1],v); p2=gtoser((GEN)x[2],v); tetpil=avma; return gerepile(av,tetpil,gdiv(p1,p2)); case t_QFR: case t_QFI: case t_VEC: case t_COL: case t_MAT: lx=lg(x); i=1; while (i= lgef(x)) return gzero; if (tx==t_SER && !signe(x)) return gzero; if (!is_recursive_t(tx)) err(talker, "this object doesn't have components (is a leaf)"); l=lontyp[tx]+n-1; if (n<1 || l>=lg(x)) err(talker,"nonexistent component"); return gcopy((GEN)x[l]); } /* with respect to the main variable if v<0, with respect to the variable v otherwise. v ignored if x is not a polynomial/series. */ GEN polcoeff0(GEN x, long n, long v) { long tx=typ(x),lx,ex,w,av,tetpil; GEN xinit; if (is_scalar_t(tx)) return n? gzero: gcopy(x); switch(tx) { case t_QFR: case t_QFI: case t_VEC: case t_COL: case t_MAT: if (n<1 || n>=lg(x)) break; return gcopy((GEN)x[n]); case t_POL: if (n<0) return gzero; w=varn(x); if (v < 0 || v == w) return (n>=lgef(x)-2)? gzero: gcopy((GEN)x[n+2]); if (v < w) return n? gzero: gcopy(x); av=avma; xinit=x; x=gsubst(gsubst(x,w,polx[MAXVARN]),v,polx[0]); if (gvar(x)) { avma=av; return n? gzero: gcopy(xinit); } if (typ(x) == t_POL) { if (n>=lgef(x)-2) { avma=av; return gzero; } x = (GEN)x[n+2]; } else x = polcoeff0(x, n, 0); tetpil=avma; return gerepile(av,tetpil,gsubst(x,MAXVARN,polx[w])); case t_SER: w=varn(x); if (v < 0 || v == w) { if (!signe(x)) return gzero; lx=lg(x); ex=valp(x); if (n=ex+lx-2) break; return gcopy((GEN)x[n-ex+2]); } if (v < w) return n? gzero: gcopy(x); av=avma; xinit=x; x=gsubst(gsubst(x,w,polx[MAXVARN]),v,polx[0]); if (gvar(x)) { avma=av; return n? gzero: gcopy(xinit); } if (gcmp0(x)) { avma=av; return gzero; } if (typ(x) == t_SER) { lx=lg(x); ex=valp(x); if (n=ex+lx-2) break; x = (GEN)x[n-ex+2]; } else x = polcoeff0(x, n, 0); tetpil=avma; return gerepile(av,tetpil,gsubst(x,MAXVARN,polx[w])); case t_RFRAC: case t_RFRACN: w = precdl; av = avma; if (v<0) v = gvar(x); ex = ggval((GEN)x[2], polx[v]); precdl = n + ex + 1; x = gtoser(x, v); precdl = w; return gerepileupto(av, polcoeff0(x, n, v)); } err(talker,"nonexistent component in truecoeff"); return NULL; /* not reached */ } GEN truecoeff(GEN x, long n) { return polcoeff0(x,n,-1); } GEN denom(GEN x) { long lx,i,av,tetpil; GEN s,t; switch(typ(x)) { case t_INT: case t_REAL: case t_INTMOD: case t_PADIC: case t_SER: return gun; case t_FRAC: case t_FRACN: return absi((GEN)x[2]); case t_COMPLEX: av=avma; t=denom((GEN)x[1]); s=denom((GEN)x[2]); tetpil=avma; return gerepile(av,tetpil,glcm(s,t)); case t_QUAD: av=avma; t=denom((GEN)x[2]); s=denom((GEN)x[3]); tetpil=avma; return gerepile(av,tetpil,glcm(s,t)); case t_POLMOD: return denom((GEN)x[2]); case t_RFRAC: case t_RFRACN: return gcopy((GEN)x[2]); case t_POL: return polun[varn(x)]; case t_VEC: case t_COL: case t_MAT: lx=lg(x); if (lx==1) return gun; av = tetpil = avma; s = denom((GEN)x[1]); for (i=2; i0)? gcopy((GEN)x[1]): gneg((GEN)x[1]); case t_POLMOD: av=avma; s=numer((GEN)x[2]); tetpil=avma; return gerepile(av,tetpil,gmodulcp(s,(GEN)x[1])); case t_RFRAC: case t_RFRACN: return gcopy((GEN)x[1]); case t_COMPLEX: case t_QUAD: case t_VEC: case t_COL: case t_MAT: av=avma; s=denom(x); tetpil=avma; return gerepile(av,tetpil,gmul(s,x)); } err(typeer,"numer"); return NULL; /* not reached */ } /* Lift only intmods if v does not occur in x, lift with respect to main * variable of x if v < 0, with respect to variable v otherwise. */ GEN lift0(GEN x, long v) { long lx,tx=typ(x),i; GEN y; switch(tx) { case t_INT: case t_REAL: return gcopy(x); case t_INTMOD: return gcopy((GEN)x[2]); case t_POLMOD: if (v < 0 || v == varn((GEN)x[1])) return gcopy((GEN)x[2]); y = cgetg(3,tx); y[1] = (long)lift0((GEN)x[1],v); y[2] = (long)lift0((GEN)x[2],v); return y; case t_SER: if (!signe(x)) return gcopy(x); lx=lg(x); y=cgetg(lx,tx); y[1]=x[1]; for (i=2; i=lontyp[tx]; i--) x[i] = (long) lift_intern0((GEN)x[i],v); return x; } err(typeer,"lift_intern"); return NULL; /* not reached */ } /* memes conventions pour v que lift */ GEN centerlift0(GEN x, long v) { long lx,tx=typ(x),i,av; GEN y; switch(tx) { case t_INT: return icopy(x); case t_INTMOD: av=avma; i=cmpii(shifti((GEN)x[2],1),(GEN)x[1]); avma=av; return (i>0)? subii((GEN)x[2],(GEN)x[1]): icopy((GEN)x[2]); case t_POLMOD: if (v < 0 || v == varn((GEN)x[1])) return gcopy((GEN)x[2]); y=cgetg(3,tx); y[1]=(long)centerlift0((GEN)x[1],v); y[2]=(long)centerlift0((GEN)x[2],v); return y; case t_SER: if (!signe(x)) return gcopy(x); lx=lg(x); y=cgetg(lx,tx); y[1]=x[1]; for (i=2; i=2; i--) if (!gcmp0(f((GEN)x[i]))) break; avma=av; if (i==1) return zeropol(varn(x)); z=cgetg(i+1,t_POL); z[1]=evalsigne(1)|evallgef(1+i)|evalvarn(varn(x)); for (j=2; j<=i; j++) z[j] = (long)f((GEN)x[j]); return z; case t_SER: if (gcmp0(x)) { z=cgetg(2,t_SER); z[1]=x[1]; return z; } lx=lg(x); av=avma; for (i=2; i=0? gun: gzero; } GEN ggt(GEN x, GEN y) { return gcmp(x,y)>0? gun: gzero; } GEN geq(GEN x, GEN y) { return _egal(x,y)? gun: gzero; } GEN gne(GEN x, GEN y) { return _egal(x,y)? gzero: gun; } GEN gand(GEN x, GEN y) { return gcmp0(x)? gzero: (gcmp0(y)? gzero: gun); } GEN gor(GEN x, GEN y) { return gcmp0(x)? (gcmp0(y)? gzero: gun): gun; } GEN gnot(GEN x) { return gcmp0(x)? gun: gzero; } /*******************************************************************/ /* */ /* FORMAL SIMPLIFICATIONS */ /* */ /*******************************************************************/ GEN geval(GEN x) { long av,tetpil,lx,i, tx = typ(x); GEN y,z; if (is_const_t(tx)) return gcopy(x); if (is_graphicvec_t(tx) || tx == t_RFRACN) { lx=lg(x); y=cgetg(lx, tx); for (i=1; ivalue; if (gegal(x, initial_value(ep))) return gcopy(z); #undef initial_value } y=gzero; av=avma; for (i=lx-1; i>1; i--) y = gadd(geval((GEN)x[i]), gmul(z,y)); return gerepileupto(av, y); case t_SER: err(impl, "evaluation of a power series"); case t_RFRAC: return gdiv(geval((GEN)x[1]),geval((GEN)x[2])); } err(typeer,"geval"); return NULL; /* not reached */ } GEN simplify_i(GEN x) { long tx=typ(x),i,lx; GEN y; switch(tx) { case t_INT: case t_REAL: case t_FRAC: case t_INTMOD: case t_PADIC: case t_QFR: case t_QFI: case t_LIST: case t_STR: case t_VECSMALL: return x; case t_FRACN: return gred(x); case t_COMPLEX: if (isexactzero((GEN)x[2])) return simplify_i((GEN)x[1]); y=cgetg(3,t_COMPLEX); y[1]=(long)simplify_i((GEN)x[1]); y[2]=(long)simplify_i((GEN)x[2]); return y; case t_QUAD: if (isexactzero((GEN)x[3])) return simplify_i((GEN)x[2]); y=cgetg(4,t_QUAD); y[1]=x[1]; y[2]=(long)simplify_i((GEN)x[2]); y[3]=(long)simplify_i((GEN)x[3]); return y; case t_POLMOD: y=cgetg(3,t_POLMOD); y[1]=(long)simplify_i((GEN)x[1]); i = typ(y[1]); if (i != t_POL) { if (i == t_INT) settyp(y, t_INTMOD); else y[1] = x[1]; /* invalid object otherwise */ } y[2]=lmod(simplify_i((GEN)x[2]),(GEN)y[1]); return y; case t_POL: lx=lgef(x); if (lx==2) return gzero; if (lx==3) return simplify_i((GEN)x[2]); y=cgetg(lx,t_POL); y[1]=x[1]; for (i=2; i 1, and scalar entries */ GEN mulmat_real(GEN x, GEN y) { long av,i,j,k, lx = lg(x), ly = lg(y), l = lg(x[1]); GEN p1, z = cgetg(ly,t_MAT); for (j=1; j=imin; i--) p1 = gadd(gmul(p1,y),(GEN)x[i]); #endif /* specific attention to sparse polynomials */ for ( ; i>=imin; i=j-1) { for (j=i; gcmp0((GEN)x[j]); j--) if (j==imin) { if (i!=j) y = gpuigs(y,i-j+1); tetpil=avma; return gerepile(av,tetpil,gmul(p1,y)); } r = (i==j)? y: gpuigs(y,i-j+1); p1 = gadd(gmul(p1,r), (GEN)x[j]); } return gerepileupto(av,p1); } p2=(GEN)x[i]; i--; r=gtrace(y); s=gneg_i(gnorm(y)); for ( ; i>=imin; i--) { p3=gadd(p2,gmul(r,p1)); p2=gadd((GEN)x[i],gmul(s,p1)); p1=p3; } p1=gmul(y,p1); tetpil=avma; return gerepile(av,tetpil,gadd(p1,p2)); }