[BACK]Return to Fgfs.c CVS log [TXT][DIR] Up to [local] / OpenXM_contrib2 / asir2000 / engine

Annotation of OpenXM_contrib2/asir2000/engine/Fgfs.c, Revision 1.16

1.16    ! noro        1: /* $OpenXM: OpenXM_contrib2/asir2000/engine/Fgfs.c,v 1.15 2003/01/04 09:06:17 noro Exp $ */
1.1       noro        2:
                      3: #include "ca.h"
                      4:
1.3       noro        5: void cont_pp_mv_sf(VL vl,VL rvl,P p,P *c,P *pp);
                      6: void gcdsf_main(VL vl,P *pa,int m,P *r);
                      7: void ugcdsf(P *pa,int m,P *r);
1.10      noro        8: void head_monomial(VL vl,V v,P p,P *coef,P *term);
1.4       noro        9: void sqfrsfmain(VL vl,P f,DCP *dcp);
                     10: void pthrootsf(P f,Q m,P *r);
                     11: void partial_sqfrsf(VL vl,V v,P f,P *r,DCP *dcp);
                     12: void gcdsf(VL vl,P *pa,int k,P *r);
1.5       noro       13: void mfctrsfmain(VL vl, P f, DCP *dcp);
1.7       noro       14: void next_evaluation_point(int *mev,int n);
1.12      noro       15: void estimatelc_sf(VL vl,VL rvl,P c,DCP dc,int *mev,P *lcp);
                     16: void mfctrsf_hensel(VL vl,VL rvl,P f,P pp0,P u0,P v0,P lcu,P lcv,int *mev,P *up);
1.7       noro       17: void substvp_sf(VL vl,VL rvl,P f,int *mev,P *r);
                     18: void shift_sf(VL vl, VL rvl, P f, int *mev, int sgn, P *r);
1.12      noro       19: void adjust_coef_sf(VL vl,VL rvl,P lcu,P u0,int *mev,P *r);
1.8       noro       20: void extended_gcd_modyk(P u0,P v0,V x,V y,int dy,P *cu,P *cv);
1.7       noro       21: void poly_to_gfsn_poly(VL vl,P f,V v,P *r);
                     22: void gfsn_poly_to_poly(VL vl,P f,V v,P *r);
1.8       noro       23: void poly_to_gfsn_poly_main(P f,V v,P *r);
                     24: void gfsn_poly_to_poly_main(P f,V v,P *r);
                     25: void gfsn_univariate_to_sfbm(P f,int dy,BM *r);
                     26: void sfbm_to_gfsn_univariate(BM f,V x,V y,P *r);
1.4       noro       27:
                     28: void lex_lc(P f,P *c)
                     29: {
                     30:        if ( !f || NUM(f) )
                     31:                *c = f;
                     32:        else
                     33:                lex_lc(COEF(DC(f)),c);
                     34: }
                     35:
                     36: DCP append_dc(DCP dc,DCP dct)
                     37: {
                     38:        DCP dcs;
                     39:
                     40:        if ( !dc )
                     41:                return dct;
                     42:        else {
                     43:                for ( dcs = dc; NEXT(dcs); dcs = NEXT(dcs) );
                     44:                NEXT (dcs) = dct;
                     45:                return dc;
                     46:        }
                     47: }
                     48:
                     49: void sqfrsf(VL vl, P f, DCP *dcp)
                     50: {
                     51:        DCP dc,dct;
                     52:        Obj obj;
1.14      noro       53:        P t,s,c,cont;
1.13      noro       54:        VL tvl,onevl;
1.4       noro       55:
                     56:        simp_ff((Obj)f,&obj); f = (P)obj;
                     57:        lex_lc(f,&c); divsp(vl,f,c,&t); f = t;
                     58:        monomialfctr(vl,f,&t,&dc); f = t;
                     59:        clctv(vl,f,&tvl); vl = tvl;
1.13      noro       60:        NEWVL(onevl); NEXT(onevl)=0;
1.4       noro       61:        if ( !vl )
                     62:                ;
                     63:        else if ( !NEXT(vl) ) {
                     64:                sfusqfr(f,&dct);
                     65:                dc = append_dc(dc,NEXT(dct));
                     66:        } else {
                     67:                t = f;
                     68:                for ( tvl = vl; tvl; tvl = NEXT(tvl) ) {
1.13      noro       69:                        onevl->v = tvl->v;
1.14      noro       70:                        cont_pp_mv_sf(vl,onevl,t,&cont,&s); t = s;
                     71:                        sqfrsf(vl,cont,&dct);
1.4       noro       72:                        dc = append_dc(dc,NEXT(dct));
                     73:                }
                     74:                sqfrsfmain(vl,t,&dct);
                     75:                dc = append_dc(dc,dct);
                     76:        }
                     77:        NEWDC(dct); DEG(dct) = ONE; COEF(dct) = (P)c; NEXT(dct) = dc;
                     78:        *dcp = dct;
                     79: }
                     80:
                     81: void sqfrsfmain(VL vl,P f,DCP *dcp)
                     82: {
                     83:        VL tvl;
                     84:        DCP dc,dct,dcs;
                     85:        P t,s;
                     86:        Q m,m1;
                     87:        V v;
                     88:
                     89:        clctv(vl,f,&tvl); vl = tvl;
                     90:        dc = 0;
                     91:        t = f;
                     92:        for ( tvl = vl; tvl; tvl = NEXT(tvl) ) {
                     93:                v = tvl->v;
                     94:                partial_sqfrsf(vl,v,t,&s,&dct); t = s;
                     95:                dc = append_dc(dc,dct);
                     96:        }
                     97:        if ( !NUM(t) ) {
                     98:                STOQ(characteristic_sf(),m);
                     99:                pthrootsf(t,m,&s);
                    100:                sqfrsfmain(vl,s,&dct);
                    101:                for ( dcs = dct; dcs; dcs = NEXT(dcs) ) {
                    102:                        mulq(DEG(dcs),m,&m1); DEG(dcs) = m1;
                    103:                }
                    104:                dc = append_dc(dc,dct);
                    105:        }
                    106:        *dcp = dc;
                    107: }
                    108:
                    109: void pthrootsf(P f,Q m,P *r)
                    110: {
                    111:        DCP dc,dc0,dct;
                    112:        N qn,rn;
                    113:
                    114:        if ( NUM(f) )
                    115:                pthrootgfs(f,r);
                    116:        else {
                    117:                dc = DC(f);
                    118:                dc0 = 0;
                    119:                for ( dc0 = 0; dc; dc = NEXT(dc) ) {
                    120:                        NEXTDC(dc0,dct);
                    121:                        pthrootsf(COEF(dc),m,&COEF(dct));
                    122:                        if ( DEG(dc) ) {
                    123:                                divn(NM(DEG(dc)),NM(m),&qn,&rn);
                    124:                                if ( rn )
                    125:                                        error("pthrootsf : cannot happen");
                    126:                                NTOQ(qn,1,DEG(dct));
                    127:                        } else
                    128:                                DEG(dct) = 0;
                    129:                }
                    130:                NEXT(dct) = 0;
                    131:                MKP(VR(f),dc0,*r);
                    132:        }
                    133: }
                    134:
                    135: void partial_sqfrsf(VL vl,V v,P f,P *r,DCP *dcp)
                    136: {
                    137:        P ps[2];
                    138:        DCP dc0,dc;
                    139:        int m;
                    140:        P t,flat,flat1,g,df,q;
                    141:
                    142:        diffp(vl,f,v,&df);
                    143:        if ( !df ) {
                    144:                *dcp = 0;
                    145:                *r = f;
                    146:                return;
                    147:        }
                    148:        ps[0] = f; ps[1] = df;
                    149:        gcdsf(vl,ps,2,&g);
                    150:        divsp(vl,f,g,&flat);
                    151:        m = 0;
                    152:        t = f;
                    153:        dc0 = 0;
                    154:        while ( !NUM(flat) ) {
                    155:                while ( divtp(vl,t,flat,&q) ) {
                    156:                        t = q; m++;
                    157:                }
                    158:                ps[0] = t; ps[1] = flat;
                    159:                gcdsf(vl,ps,2,&flat1);
                    160:                divsp(vl,flat,flat1,&g);
                    161:                flat = flat1;
                    162:                NEXTDC(dc0,dc);
                    163:                COEF(dc) = g;
                    164:                STOQ(m,DEG(dc));
                    165:        }
                    166:        NEXT(dc) = 0;
                    167:        *dcp = dc0;
                    168:        *r = t;
                    169: }
1.1       noro      170:
                    171: void gcdsf(VL vl,P *pa,int k,P *r)
                    172: {
1.3       noro      173:        P *ps,*pl,*pm;
                    174:        P **cp;
1.1       noro      175:        int *cn;
                    176:        DCP *ml;
                    177:        Obj obj;
                    178:        int i,j,l,m;
                    179:        P mg,mgsf,t;
                    180:        VL avl,nvl,tvl,svl;
                    181:
                    182:        ps = (P *)ALLOCA(k*sizeof(P));
                    183:        for ( i = 0, m = 0; i < k; i++ ) {
                    184:                simp_ff((Obj)pa[i],&obj);
                    185:                if ( obj )
1.3       noro      186:                        ps[m++] = (P)obj;
1.1       noro      187:        }
                    188:        if ( !m ) {
                    189:                *r = 0;
                    190:                return;
                    191:        }
                    192:        if ( m == 1 ) {
1.3       noro      193:                *r = ps[0];
1.1       noro      194:                return;
                    195:        }
                    196:        pl = (P *)ALLOCA(m*sizeof(P));
                    197:        ml = (DCP *)ALLOCA(m*sizeof(DCP));
                    198:        for ( i = 0; i < m; i++ )
                    199:                monomialfctr(vl,ps[i],&pl[i],&ml[i]);
1.3       noro      200:        gcdmonomial(vl,ml,m,&mg); simp_ff((Obj)mg,&obj); mgsf = (P)obj;
1.1       noro      201:        for ( i = 0, nvl = vl, avl = 0; nvl && i < m; i++ ) {
                    202:                clctv(vl,pl[i],&tvl);
                    203:                intersectv(nvl,tvl,&svl); nvl = svl;
                    204:                mergev(vl,avl,tvl,&svl); avl = svl;
                    205:        }
                    206:        if ( !nvl ) {
                    207:                *r = mgsf;
                    208:                return;
                    209:        }
                    210:        if ( !NEXT(avl) ) {
                    211:                ugcdsf(pl,m,&t);
                    212:                mulp(vl,mgsf,t,r);
                    213:                return;
                    214:        }
                    215:        for ( tvl = nvl, i = 0; tvl; tvl = NEXT(tvl), i++ );
                    216:        for ( tvl = avl, j = 0; tvl; tvl = NEXT(tvl), j++ );
                    217:        if ( i == j ) {
                    218:                /* all the pl[i]'s have the same variables */
                    219:                gcdsf_main(avl,pl,m,&t);
                    220:        } else {
                    221:                cp = (P **)ALLOCA(m*sizeof(P *));
                    222:                cn = (int *)ALLOCA(m*sizeof(int));
                    223:                for ( i = 0; i < m; i++ ) {
                    224:                        cp[i] = (P *)ALLOCA(lengthp(pl[i])*sizeof(P));
                    225:                        cn[i] = pcoef(vl,nvl,pl[i],cp[i]);
                    226:                }
                    227:                for ( i = j = 0; i < m; i++ )
                    228:                        j += cn[i];
                    229:                pm = (P *)ALLOCA(j*sizeof(P));
                    230:                for ( i = l = 0; i < m; i++ )
                    231:                        for ( j = 0; j < cn[i]; j++ )
                    232:                                pm[l++] = cp[i][j];
                    233:                gcdsf(vl,pm,l,&t);
                    234:        }
                    235:        mulp(vl,mgsf,t,r);
                    236: }
                    237:
                    238: /* univariate gcd */
                    239:
                    240: void ugcdsf(P *pa,int m,P *r)
                    241: {
1.3       noro      242:        P *ps;
1.1       noro      243:        int i;
                    244:        UM w1,w2,w3,w;
                    245:        int d;
                    246:        V v;
                    247:
                    248:        if ( m == 1 ) {
                    249:                *r = pa[0];
                    250:                return;
                    251:        }
1.3       noro      252:        for ( i = 0; i < m; i++ )
                    253:                if ( NUM(pa[i]) ) {
                    254:                        itogfs(1,r);
                    255:                        return;
                    256:                }
1.1       noro      257:        ps = (P *)ALLOCA(m*sizeof(P));
                    258:        sort_by_deg(m,pa,ps);
1.3       noro      259:        v = VR(ps[m-1]);
                    260:        d = getdeg(v,ps[m-1]);
1.1       noro      261:        w1 = W_UMALLOC(d);
                    262:        w2 = W_UMALLOC(d);
                    263:        w3 = W_UMALLOC(d);
                    264:        ptosfum(ps[0],w1);
                    265:        for ( i = 1; i < m; i++ ) {
                    266:                ptosfum(ps[i],w2);
                    267:                gcdsfum(w1,w2,w3);
                    268:                w = w1; w1 = w3; w3 = w;
                    269:                if ( !DEG(w1) ) {
1.3       noro      270:                        itogfs(1,r);
1.1       noro      271:                        return;
                    272:                }
                    273:        }
                    274:        sfumtop(v,w1,r);
                    275: }
                    276:
1.4       noro      277: /* deg(HT(p),v), where p is considered as distributed poly over F[v] */
                    278: int gethdeg(VL vl,V v,P p)
                    279: {
                    280:        DCP dc;
                    281:        Q dmax;
                    282:        P cmax;
                    283:
                    284:        if ( !p )
                    285:                return -1;
                    286:        else if ( NUM(p) )
                    287:                return 0;
                    288:        else if ( VR(p) != v )
                    289:                /* HT(p) = HT(lc(p))*x^D */
                    290:                return gethdeg(vl,v,COEF(DC(p)));
                    291:        else {
                    292:                /* VR(p) = v */
                    293:                dc = DC(p); dmax = DEG(dc); cmax = COEF(dc);
                    294:                for ( dc = NEXT(dc); dc; dc = NEXT(dc) )
                    295:                        if ( compp(vl,COEF(dc),cmax) > 0 ) {
                    296:                                dmax = DEG(dc); cmax = COEF(dc);
                    297:                        }
                    298:                return QTOS(dmax);
                    299:        }
                    300: }
1.1       noro      301:
                    302: /* all the pa[i]'s have the same variables (=vl) */
                    303:
                    304: void gcdsf_main(VL vl,P *pa,int m,P *r)
                    305: {
1.3       noro      306:        int nv,i,i0,imin,d,d0,d1,d2,dmin,index;
                    307:        V v,v0,vmin;
1.2       noro      308:        VL tvl,nvl,rvl,nvl0,rvl0;
1.3       noro      309:        P *pc, *ps, *ph,*lps;
                    310:        P x,t,cont,hg,g,hm,mod,s;
                    311:        P hge,ge,ce,he,u,cof1e,mode,mod1,adj,cof1,coadj,q;
                    312:        GFS sf;
1.2       noro      313:
1.1       noro      314:        for ( nv = 0, tvl = vl; tvl; tvl = NEXT(tvl), nv++);
                    315:        if ( nv == 1 ) {
                    316:                ugcdsf(pa,m,r);
                    317:                return;
                    318:        }
1.4       noro      319:        /* find v s.t. min(deg(pa[i],v)+gethdeg(pa[i],v)) is minimal */
1.1       noro      320:        tvl = vl;
                    321:        do {
                    322:                v = tvl->v;
                    323:                i = 0;
                    324:                do {
1.4       noro      325:                        d = getdeg(v,pa[i])+gethdeg(vl,v,pa[i]);
1.1       noro      326:                        if ( i == 0 || (d < d0) ) {
                    327:                                d0 = d; i0 = i; v0 = v;
                    328:                        }
                    329:                } while ( ++i < m );
                    330:                if ( tvl == vl || (d0 < dmin) ) {
                    331:                        dmin = d0; imin = i0; vmin = v0;
                    332:                }
                    333:        } while ( tvl = NEXT(tvl) );
                    334:
                    335:        /* reorder variables so that vmin is the last variable */
                    336:        for ( nvl0 = 0, rvl0 = 0, tvl = vl; tvl; tvl = NEXT(tvl) )
                    337:                if ( tvl->v != vmin ) {
                    338:                        NEXTVL(nvl0,nvl); nvl->v = tvl->v;
                    339:                        NEXTVL(rvl0,rvl); rvl->v = tvl->v;
                    340:                }
                    341:        /* rvl = remaining variables */
1.3       noro      342:        NEXT(rvl) = 0; rvl = rvl0;
1.1       noro      343:        /* nvl = ...,vmin */
1.3       noro      344:        NEXTVL(nvl0,nvl); nvl->v = vmin; NEXT(nvl) = 0; nvl = nvl0;
1.2       noro      345:        MKV(vmin,x);
1.1       noro      346:
                    347:        /* for content and primitive part */
                    348:        pc = (P *)ALLOCA(m*sizeof(P));
                    349:        ps = (P *)ALLOCA(m*sizeof(P));
                    350:        ph = (P *)ALLOCA(m*sizeof(P));
                    351:        /* separate the contents */
                    352:        for ( i = 0; i < m; i++ ) {
1.3       noro      353:                reorderp(nvl,vl,pa[i],&t);
1.1       noro      354:                cont_pp_mv_sf(nvl,rvl,t,&pc[i],&ps[i]);
1.10      noro      355:                head_monomial(nvl,vmin,ps[i],&ph[i],&t);
1.1       noro      356:        }
                    357:        ugcdsf(pc,m,&cont);
                    358:        ugcdsf(ph,m,&hg);
                    359:
                    360:        /* for hg*pp (used in check phase) */
                    361:        lps = (P *)ALLOCA(m*sizeof(P));
                    362:        for ( i = 0; i < m; i++ )
                    363:                mulp(nvl,hg,ps[i],&lps[i]);
                    364:
                    365:        while ( 1 ) {
                    366:                g = 0;
1.3       noro      367:                cof1 = 0;
1.1       noro      368:                hm = 0;
1.3       noro      369:                itogfs(1,&mod);
1.1       noro      370:                index = 0;
1.3       noro      371:                for ( index = 0; getdeg(vmin,mod) <= d+1; index++ ) {
1.1       noro      372:                        /* evaluation pt */
1.3       noro      373:                        indextogfs(index,&s);
1.1       noro      374:                        substp(nvl,hg,vmin,s,&hge);
                    375:                        if ( !hge )
                    376:                                continue;
                    377:                        for ( i = 0; i < m; i++ )
                    378:                                substp(nvl,ps[i],vmin,s,&ph[i]);
                    379:                        /* ge = GCD(ps[0]|x=s,...,ps[m-1]|x=s) */
                    380:                        gcdsf(nvl,ph,m,&ge);
1.10      noro      381:                        head_monomial(nvl,vmin,ge,&ce,&he);
1.3       noro      382:                        if ( NUM(he) ) {
1.1       noro      383:                                *r = cont;
                    384:                                return;
                    385:                        }
1.3       noro      386:                        divgfs((GFS)hge,(GFS)ce,&sf); t = (P)sf;
                    387:                        mulp(nvl,t,ge,&u); ge = u;
1.1       noro      388:                        divsp(nvl,ph[imin],ge,&t); mulp(nvl,hge,t,&cof1e);
1.2       noro      389:                        /* hm=0 : reset; he==hm : lucky */
1.3       noro      390:                        if ( !hm || !compp(nvl,he,hm) ) {
1.2       noro      391:                                substp(nvl,mod,vmin,s,&mode); divsp(nvl,mod,mode,&mod1);
                    392:                                /* adj = mod/(mod|x=s)*(ge-g|x=s) */
                    393:                                substp(nvl,g,vmin,s,&t);
                    394:                                subp(nvl,ge,t,&u); mulp(nvl,mod1,u,&adj);
                    395:                                /* coadj = mod/(mod|vmin=s)*(cof1e-cof1e|vmin=s) */
                    396:                                substp(nvl,cof1,vmin,s,&t);
1.3       noro      397:                                subp(nvl,cof1e,t,&u); mulp(nvl,mod1,u,&coadj);
1.2       noro      398:                                if ( !adj ) {
                    399:                                        /* adj == gcd ? */
                    400:                                        for ( i = 0; i < m; i++ )
1.3       noro      401:                                                if ( !divtp(nvl,lps[i],g,&t) )
1.2       noro      402:                                                        break;
                    403:                                        if ( i == m ) {
1.3       noro      404:                                                cont_pp_mv_sf(nvl,rvl,g,&t,&u);
1.2       noro      405:                                                mulp(nvl,cont,u,&t);
1.3       noro      406:                                                reorderp(vl,nvl,t,r);
1.2       noro      407:                                                return;
                    408:                                        }
                    409:                                } else if ( !coadj ) {
1.3       noro      410:                                        /* ps[imin]/coadj == gcd ? */
                    411:                                        if ( divtp(nvl,lps[imin],cof1,&q) ) {
1.2       noro      412:                                                for ( i = 0; i < m; i++ )
                    413:                                                        if ( !divtp(nvl,lps[i],q,&t) )
                    414:                                                                break;
                    415:                                                if ( i == m ) {
                    416:                                                        cont_pp_mv_sf(nvl,rvl,q,&t,&u);
                    417:                                                        mulp(nvl,cont,u,&t);
1.3       noro      418:                                                        reorderp(vl,nvl,t,r);
1.2       noro      419:                                                        return;
                    420:                                                }
                    421:                                        }
                    422:                                }
                    423:                                addp(nvl,g,adj,&t); g = t;
                    424:                                addp(nvl,cof1,coadj,&t); cof1 = t;
                    425:                                subp(nvl,x,s,&t); mulp(nvl,mod,t,&u); mod = u;
                    426:                                hm = he;
                    427:                        } else {
                    428:                                d1 = homdeg(hm); d2 = homdeg(he);
                    429:                                if ( d1 < d2 ) /* we use current hm */
                    430:                                        continue;
                    431:                                else if ( d1 > d2 ) {
                    432:                                        /* use he */
                    433:                                        g = ge;
                    434:                                        cof1 = cof1e;
                    435:                                        hm = he;
                    436:                                        subp(nvl,x,s,&mod);
                    437:                                } else {
                    438:                                        /* d1==d2, but hm!=he => both are unlucky */
                    439:                                        g = 0;
                    440:                                        cof1 = 0;
1.3       noro      441:                                        itogfs(1,&mod);
1.2       noro      442:                                }
1.1       noro      443:                        }
                    444:                }
                    445:        }
                    446: }
                    447:
1.10      noro      448: void head_monomial(VL vl,V v,P p,P *coef,P *term)
1.1       noro      449: {
                    450:        P t,s,u;
                    451:        DCP dc;
                    452:        GFS one;
                    453:
1.3       noro      454:        itogfs(1,&one);
                    455:        t = (P)one;
1.1       noro      456:        while ( 1 ) {
                    457:                if ( NUM(p) || VR(p) == v ) {
                    458:                        *coef = p;
                    459:                        *term = t;
                    460:                        return;
                    461:                } else {
1.3       noro      462:                        NEWDC(dc);
                    463:                        COEF(dc) = (P)one; DEG(dc) = DEG(DC(p));
1.1       noro      464:                        MKP(VR(p),dc,s);
                    465:                        mulp(vl,t,s,&u); t = u;
                    466:                        p = COEF(DC(p));
                    467:                }
                    468:        }
                    469: }
                    470:
                    471: void cont_pp_mv_sf(VL vl,VL rvl,P p,P *c,P *pp)
                    472: {
                    473:        DP dp;
                    474:        MP t;
                    475:        int i,m;
                    476:        P *ps;
1.16    ! noro      477:        struct order_spec spec, currentspec;
        !           478:        extern struct order_spec dp_current_spec;
1.1       noro      479:
1.16    ! noro      480:        currentspec = dp_current_spec;
1.15      noro      481:        create_order_spec(0,&spec);
                    482:        initd(&spec);
1.1       noro      483:        ptod(vl,rvl,p,&dp);
                    484:        for ( t = BDY(dp), m = 0; t; t = NEXT(t), m++ );
                    485:        ps = (P *)ALLOCA(m*sizeof(P));
1.3       noro      486:        for ( t = BDY(dp), i = 0; t; t = NEXT(t), i++ )
1.1       noro      487:                ps[i] = C(t);
1.10      noro      488:        gcdsf(vl,ps,m,c);
1.3       noro      489:        divsp(vl,p,*c,pp);
1.16    ! noro      490:        initd(&currentspec);
1.5       noro      491: }
                    492:
                    493: void mfctrsf(VL vl, P f, DCP *dcp)
                    494: {
                    495:        DCP dc0,dc,dct,dcs,dcr;
                    496:        Obj obj;
                    497:
                    498:        simp_ff((Obj)f,&obj); f = (P)obj;
                    499:        sqfrsf(vl,f,&dct);
                    500:        dc = dc0 = dct; dct = NEXT(dct); NEXT(dc) = 0;
                    501:        for ( ; dct; dct = NEXT(dct) ) {
                    502:                mfctrsfmain(vl,COEF(dct),&dcs);
                    503:                for ( dcr = dcs; dcr; dcr = NEXT(dcr) )
                    504:                        DEG(dcr) = DEG(dct);
                    505:                for ( ; NEXT(dc); dc = NEXT(dc) );
                    506:                NEXT(dc) = dcs;
                    507:        }
                    508:        *dcp = dc0;
                    509: }
                    510:
                    511: /* f : sqfr, non const */
                    512:
                    513: void mfctrsfmain(VL vl, P f, DCP *dcp)
                    514: {
1.6       noro      515:        VL tvl,nvl,rvl;
1.7       noro      516:        DCP dc,dc0,dc1,dc2,dct,lcfdc,dcs;
                    517:        int imin,inext,i,j,n,k,np;
1.5       noro      518:        int *da;
                    519:        V vx,vy;
                    520:        V *va;
1.7       noro      521:        P *l,*tl;
1.5       noro      522:        P gcd,g,df,dfmin;
                    523:        P pa[2];
1.10      noro      524:        P f0,pp0,spp0,c,c0,x,y,u,v,lcf,lcu,lcv,u0,v0,t,s;
1.12      noro      525:        P ype,yme,fin;
1.6       noro      526:        GFS ev,evy;
                    527:        P *fp0;
                    528:        int *mev,*win;
1.5       noro      529:
                    530:        clctv(vl,f,&tvl); vl = tvl;
                    531:        if ( !vl )
                    532:                error("mfctrsfmain : cannot happen");
                    533:        if ( !NEXT(vl) ) {
                    534:                /* univariate */
                    535:                ufctrsf(f,&dc);
                    536:                /* remove lc */
                    537:                *dcp = NEXT(dc);
                    538:                return;
                    539:        }
                    540:        for ( n = 0, tvl = vl; tvl; tvl = NEXT(tvl), n++ );
                    541:        va = (V *)ALLOCA(n*sizeof(int));
                    542:        da = (int *)ALLOCA(n*sizeof(int));
                    543:        /* find v s.t. diff(f,v) is nonzero and deg(f,v) is minimal */
                    544:        imin = -1;
                    545:        for ( i = 0, tvl = vl; i < n; tvl = NEXT(tvl), i++ ) {
                    546:                va[i] = tvl->v;
                    547:                da[i] = getdeg(va[i],f);
                    548:                diffp(vl,f,va[i],&df);
                    549:                if ( !df )
                    550:                        continue;
                    551:                if ( imin < 0 || da[i] < da[imin] ) {
                    552:                        dfmin = df;
                    553:                        imin = i;
                    554:                }
                    555:        }
                    556:        /* find v1 neq v s.t. deg(f,v) is minimal */
                    557:        inext = -1;
                    558:        for ( i = 0; i < n; i++ ) {
                    559:                if ( i == imin )
                    560:                        continue;
                    561:                if ( inext < 0 || da[i] < da[inext] )
                    562:                        inext = i;
                    563:        }
                    564:        pa[0] = f;
                    565:        pa[1] = dfmin;
1.11      noro      566:        gcdsf(vl,pa,2,&gcd);
1.5       noro      567:        if ( !NUM(gcd) ) {
                    568:                /* f = gcd * f/gcd */
                    569:                mfctrsfmain(vl,gcd,&dc1);
                    570:                divsp(vl,f,gcd,&g);
                    571:                mfctrsfmain(vl,g,&dc2);
                    572:                for ( dct = dc1; NEXT(dct); dct = NEXT(dct) );
                    573:                NEXT(dct) = dc2;
                    574:                *dcp = dc1;
                    575:                return;
                    576:        }
                    577:        /* create vl s.t. vl[0] = va[imin], vl[1] = va[inext] */
                    578:        nvl = 0;
                    579:        NEXTVL(nvl,tvl); tvl->v = va[imin];
                    580:        NEXTVL(nvl,tvl); tvl->v = va[inext];
                    581:        for ( i = 0; i < n; i++ ) {
                    582:                if ( i == imin || i == inext )
                    583:                        continue;
                    584:                NEXTVL(nvl,tvl); tvl->v = va[i];
                    585:        }
                    586:        NEXT(tvl) = 0;
                    587:
1.12      noro      588:        fin = f;
1.10      noro      589:        reorderp(nvl,vl,f,&g); f = g;
1.5       noro      590:        vx = nvl->v;
                    591:        vy = NEXT(nvl)->v;
1.6       noro      592:        MKV(vx,x);
                    593:        MKV(vy,y);
                    594:        /* remaining variables */
                    595:        rvl = NEXT(NEXT(nvl));
                    596:        if ( !rvl ) {
1.5       noro      597:                /* bivariate */
1.10      noro      598:                sfbfctr(f,vx,vy,getdeg(vx,f),&dc1);
1.5       noro      599:                for ( dc0 = 0; dc1; dc1 = NEXT(dc1) ) {
                    600:                        NEXTDC(dc0,dc);
                    601:                        DEG(dc) = ONE;
                    602:                        reorderp(vl,nvl,COEF(dc1),&COEF(dc));
                    603:                }
                    604:                NEXT(dc) = 0;
                    605:                *dcp = dc0;
                    606:                return;
                    607:        }
1.6       noro      608:        /* n >= 3;  nvl = (vx,vy,X) */
                    609:        /* find good evaluation pt for X */
                    610:        mev = (int *)CALLOC(n-2,sizeof(int));
                    611:        while ( 1 ) {
1.10      noro      612:                /* lcf(mev)=0 => invalid */
                    613:                substvp_sf(nvl,rvl,COEF(DC(f)),mev,&t);
                    614:                if ( t ) {
                    615:                        substvp_sf(nvl,rvl,f,mev,&f0);
                    616:                        pa[0] = f0;
                    617:                        diffp(nvl,f0,vx,&pa[1]);
                    618:                        if ( pa[1] ) {
                    619:                                gcdsf(nvl,pa,2,&gcd);
1.6       noro      620:                        /* XXX maybe we have to accept the case where gcd is a poly of y */
1.10      noro      621:                                if ( NUM(gcd) )
                    622:                                        break;
                    623:                        }
1.6       noro      624:                }
1.7       noro      625:                /* XXX if generated indices exceed q of GF(q) => error in indextogfs */
                    626:                next_evaluation_point(mev,n-2);
1.6       noro      627:        }
1.10      noro      628:        /* f0 = f(x,y,mev) */
                    629:        /* separate content; f0 may have the content wrt x */
                    630:        cont_pp_sfp(nvl,f0,&c0,&pp0);
1.6       noro      631:
1.7       noro      632:        /* factorize pp0; pp0 = pp0(x,y+evy) = prod dc */
                    633:        sfbfctr_shift(pp0,vx,vy,getdeg(vx,pp0),&evy,&spp0,&dc); pp0 = spp0;
1.6       noro      634:
                    635:        if ( !NEXT(dc) ) {
                    636:                /* f is irreducible */
1.12      noro      637:                NEWDC(dc); DEG(dc) = ONE; COEF(dc) = fin; NEXT(dc) = 0;
1.6       noro      638:                *dcp = dc;
                    639:                return;
                    640:        }
1.7       noro      641:        /* ype = y+evy, yme = y-evy */
                    642:        addp(nvl,y,(P)evy,&ype); subp(nvl,y,(P)evy,&yme);
                    643:
1.6       noro      644:        /* shift c0; c0 <- c0(y+evy) */
1.7       noro      645:        substp(nvl,c0,vy,ype,&s); c0 = s;
                    646:
                    647:        /* shift f; f <- f(y+evy) */
                    648:        substp(nvl,f,vy,ype,&s); f = s;
                    649:
                    650:        /* now f(x,0,mev) = c0 * prod dc */
1.6       noro      651:
                    652:        /* factorize lc_x(f) */
                    653:        lcf = COEF(DC(f));
1.7       noro      654:        mfctrsf(nvl,lcf,&dct);
                    655:        /* skip the first element (= a number) */
1.12      noro      656:        lcfdc = NEXT(dct);
1.6       noro      657:
                    658:        /* np = number of bivariate factors */
                    659:        for ( np = 0, dct = dc; dct; dct = NEXT(dct), np++ );
                    660:        fp0 = (P *)ALLOCA((np+1)*sizeof(P));
                    661:        for ( i = 0, dct = dc; i < np; dct = NEXT(dct), i++ )
                    662:                fp0[i] = COEF(dct);
                    663:        fp0[np] = 0;
1.7       noro      664:        l = tl = (P *)ALLOCA((np+1)*sizeof(P));
1.6       noro      665:        win = W_ALLOC(np+1);
1.7       noro      666:
1.6       noro      667:        for ( k = 1, win[0] = 1, --np; ; ) {
                    668:                itogfs(1,&u0);
                    669:                /* u0 = product of selected factors */
                    670:                for ( i = 0; i < k; i++ ) {
                    671:                        mulp(nvl,u0,fp0[win[i]],&t); u0 = t;
                    672:                }
                    673:                /* we have to consider the content */
1.10      noro      674:                /* f0 = c0*u0*v0 */
1.12      noro      675:                mulp(nvl,LC(u0),c0,&c); estimatelc_sf(nvl,rvl,c,lcfdc,mev,&lcu);
1.6       noro      676:                divsp(nvl,pp0,u0,&v0);
1.12      noro      677:                mulp(nvl,LC(v0),c0,&c); estimatelc_sf(nvl,rvl,c,lcfdc,mev,&lcv);
                    678:                mfctrsf_hensel(nvl,rvl,f,pp0,u0,v0,lcu,lcv,mev,&u);
1.7       noro      679:                if ( u ) {
                    680:                        /* save the factor */
                    681:                        reorderp(vl,nvl,u,&t);
1.12      noro      682:                        /* y -> y-evy */
                    683:                        substp(vl,t,vy,yme,tl++);
1.7       noro      684:
                    685:                        /* update f,pp0 */
                    686:                        divsp(nvl,f,u,&t); f = t;
                    687:                        divsp(nvl,pp0,u0,&t); pp0 = t;
                    688:                        /* update win, fp0 */
                    689:                        for ( i = 0; i < k-1; i++ )
                    690:                        for ( j = win[i]+1; j < win[i+1]; j++ )
                    691:                                fp0[j-i-1] = fp0[j];
                    692:                        for ( j = win[k-1]+1; j <= np; j++ )
                    693:                                        fp0[j-k] = fp0[j];
                    694:                        if ( ( np -= k ) < k ) break;
                    695:                        if ( np-win[0]+1 < k )
                    696:                                if ( ++k <= np ) {
                    697:                                        for ( i = 0; i < k; i++ )
                    698:                                                win[i] = i + 1;
                    699:                                        continue;
                    700:                                } else
                    701:                                        break;
                    702:                        else
                    703:                                for ( i = 1; i < k; i++ )
                    704:                                        win[i] = win[0] + i;
                    705:                } else {
                    706:                        if ( ncombi(1,np,k,win) == 0 )
                    707:                                if ( k == np ) break;
                    708:                                else
                    709:                                        for ( i = 0, ++k; i < k; i++ )
                    710:                                                win[i] = i + 1;
                    711:                }
1.6       noro      712:        }
1.10      noro      713:        reorderp(vl,nvl,f,&t);
1.12      noro      714:        /* y -> y-evy */
                    715:        substp(vl,t,vy,yme,tl++);
1.10      noro      716:        *tl = 0;
                    717:        for ( dc0 = 0, i = 0; l[i]; i++ ) {
                    718:                NEXTDC(dc0,dc); DEG(dc) = ONE; COEF(dc) = l[i];
                    719:        }
                    720:        NEXT(dc) = 0; *dcp = dc0;
1.6       noro      721: }
                    722:
1.7       noro      723: void next_evaluation_point(int *e,int n)
                    724: {
                    725:        int i,t,j;
                    726:
                    727:        for ( i = n-1; i >= 0; i-- )
                    728:                if ( e[i] ) break;
                    729:        if ( i < 0 ) e[n-1] = 1;
                    730:        else if ( i == 0 ) {
                    731:                t = e[0]; e[0] = 0; e[n-1] = t+1;
                    732:        } else {
                    733:                e[i-1]++; t = e[i];
                    734:                for ( j = i; j < n-1; j++ )
                    735:                        e[j] = 0;
                    736:                e[n-1] = t-1;
                    737:        }
                    738: }
                    739:
                    740: /*
                    741:  * dc : f1^E1*...*fk^Ek
1.12      noro      742:  * find e1,...,ek s.t. fi(mev)^ei | c
1.7       noro      743:  * and return f1^e1*...*fk^ek
                    744:  * vl = (vx,vy,rvl)
                    745:  */
                    746:
1.12      noro      747: void estimatelc_sf(VL vl,VL rvl,P c,DCP dc,int *mev,P *lcp)
1.7       noro      748: {
                    749:        DCP dct;
                    750:        P r,c1,c2,t,s,f;
                    751:        int i,d;
                    752:        Q q;
                    753:
                    754:        for ( dct = dc, r = (P)ONE; dct; dct = NEXT(dct) ) {
                    755:                if ( NUM(COEF(dct)) )
                    756:                        continue;
                    757:                /* constant part */
1.12      noro      758:                substvp_sf(vl,rvl,COEF(dct),mev,&f);
1.7       noro      759:                d = QTOS(DEG(dct));
                    760:                for ( i = 0, c1 = c; i < d; i++ )
                    761:                        if ( !divtp(vl,c1,f,&c2) )
                    762:                                break;
                    763:                        else
                    764:                                c1 = c2;
                    765:                if ( i ) {
                    766:                        STOQ(i,q);
                    767:                        pwrp(vl,COEF(dct),q,&s); mulp(vl,r,s,&t); r = t;
                    768:                }
                    769:        }
                    770:        *lcp = r;
                    771: }
                    772:
                    773: void substvp_sf(VL vl,VL rvl,P f,int *mev,P *r)
                    774: {
                    775:        int i;
                    776:        VL tvl;
                    777:        P g,t;
                    778:        GFS ev;
                    779:
                    780:        for ( g = f, i = 0, tvl = rvl; tvl; tvl = NEXT(tvl), i++ ) {
                    781:                if ( !mev )
                    782:                        ev = 0;
                    783:                else
                    784:                        indextogfs(mev[i],&ev);
                    785:                substp(vl,g,tvl->v,(P)ev,&t); g = t;
                    786:        }
                    787:        *r = g;
                    788: }
                    789:
                    790: /*
                    791:  * f <- f(X+sgn*mev)
                    792:  */
                    793:
                    794: void shift_sf(VL vl, VL rvl, P f, int *mev, int sgn, P *r)
                    795: {
                    796:        VL tvl;
                    797:        int i;
                    798:        P x,g,t,s;
                    799:        GFS ev;
                    800:
                    801:        for ( g = f, tvl = rvl, i = 0; tvl; tvl = NEXT(tvl), i++ ) {
                    802:                if ( !mev[i] )
                    803:                        continue;
                    804:                indextogfs(mev[i],&ev);
                    805:                MKV(tvl->v,x);
                    806:                if ( sgn > 0 )
                    807:                        addp(vl,x,(P)ev,&t);
                    808:                else
                    809:                        subp(vl,x,(P)ev,&t);
                    810:                substp(vl,g,tvl->v,t,&s); g = s;
                    811:        }
                    812:        *r = g;
                    813: }
                    814:
                    815: /*
                    816:  * pp(f(0)) = u0*v0
                    817:  */
                    818:
1.12      noro      819: void mfctrsf_hensel(VL vl,VL rvl,P f,P pp0,P u0,P v0,P lcu,P lcv,int *mev,P *up)
1.7       noro      820: {
                    821:        VL tvl,onevl;
                    822:        P t,s,w,u,v,ff,si,wu,wv,fj,cont;
                    823:        UM ydy;
                    824:        V vx,vy;
                    825:        int dy,n,i,dbd,nv,j;
                    826:        int *md;
                    827:        P *uh,*vh;
1.12      noro      828:        P x,du0,dv0,m,q,r,fin;
1.7       noro      829:        P *cu,*cv;
                    830:        GFSN inv;
                    831:
1.13      noro      832:        /* check the validity of lc's and adjust coeffs */
                    833:        /* f                -> lcu*lcv*x^(m+l)+... */
                    834:        mulp(vl,lcu,lcv,&t);
                    835:        if ( !divtp(vl,t,LC(f),&m) ) {
                    836:                *up = 0; return;
                    837:        }
                    838:        mulp(vl,m,f,&t); f = t;
1.12      noro      839:        /* u0 = am x^m+ ... -> lcu*x^m + a(m-1)*(lcu(mev)/am)*x^(m-1)+... */
                    840:        /* v0 = bm x^l+ ... -> lcv*x^l + b(l-1)*(lcv(mev)/bl)*x^(l-1)+... */
                    841:        adjust_coef_sf(vl,rvl,lcu,u0,mev,&u);
                    842:        adjust_coef_sf(vl,rvl,lcv,v0,mev,&v);
1.10      noro      843:
1.12      noro      844:        /* f <- f(X+mev), u <- u(X+mev), v <- v(X+mev) */
                    845:        fin = f;
                    846:        shift_sf(vl,rvl,f,mev,1,&s); f = s;
                    847:        shift_sf(vl,rvl,u,mev,1,&s); u = s;
                    848:        shift_sf(vl,rvl,v,mev,1,&s); v = s;
                    849:
1.7       noro      850:        vx = vl->v; vy = NEXT(vl)->v;
                    851:        n = getdeg(vx,f);
                    852:        dy = getdeg(vy,f)+1;
                    853:        MKV(vx,x);
                    854:        cu = (P *)ALLOCA((n+1)*sizeof(P));
                    855:        cv = (P *)ALLOCA((n+1)*sizeof(P));
                    856:
                    857:        /* ydy = y^dy */
1.10      noro      858:        ydy = C_UMALLOC(dy); DEG(ydy) = dy; COEF(ydy)[dy] = _onesf();
1.7       noro      859:        setmod_gfsn(ydy);
                    860:
                    861:        /* (R[y]/(y^dy))[x,X] */
1.10      noro      862:        poly_to_gfsn_poly(vl,f,vy,&ff);
1.7       noro      863:        poly_to_gfsn_poly(vl,u,vy,&t); u = t;
                    864:        poly_to_gfsn_poly(vl,v,vy,&t); v = t;
                    865:        substvp_sf(vl,rvl,u,0,&u0);
                    866:        substvp_sf(vl,rvl,v,0,&v0);
                    867:
                    868:        /* compute a(x,y), b(x,y) s.t. a*u0+b*v0 = 1 mod y^dy */
1.8       noro      869:        extended_gcd_modyk(u0,v0,vx,vy,dy,&cu[0],&cv[0]);
1.7       noro      870:
                    871:        /* dv0 = LC(v0)^(-1)*v0 mod y^dy */
                    872:        invgfsn((GFSN)LC(v0),&inv); mulp(vl,v0,(P)inv,&dv0);
                    873:
                    874:        /* cu[i]*u0+cv[i]*v0 = x^i mod y^dy */
1.10      noro      875:        /* (x*cu[i])*u0+(x*cv[i])*v0 = x^(i+1) */
                    876:        /* x*cu[i] = q*dv0+r => cu[i+1] = r */
                    877:        /* cv[i+1]*v0 = x*cv[i]*v0+q*u0*dv0 = (x*cv[i]+q*u0*inv)*v0 */
1.7       noro      878:        for ( i = 1; i <= n; i++ ) {
                    879:                mulp(vl,x,cu[i-1],&m); divsrp(vl,m,dv0,&q,&cu[i]);
1.10      noro      880:                mulp(vl,x,cv[i-1],&m); mulp(vl,q,(P)inv,&t);
                    881:                mulp(vl,t,u0,&s);
                    882:                addp(vl,m,s,&cv[i]);
                    883:        }
                    884:
                    885: #if 0
                    886:        /* XXX : check */
                    887:        for ( i = 0; i <= n; i++ ) {
                    888:                mulp(vl,cu[i],u0,&m); mulp(vl,cv[i],v0,&s);
                    889:                addp(vl,m,s,&w);
                    890:                printexpr(vl,w);
                    891:                fprintf(asir_out,"\n");
1.7       noro      892:        }
1.10      noro      893: #endif
                    894:
1.7       noro      895:        dbd = dbound(vx,f)+1;
                    896:
                    897:        /* extract homogeneous parts */
                    898:        W_CALLOC(dbd,P,uh); W_CALLOC(dbd,P,vh);
                    899:        for ( i = 0; i <= dbd; i++ ) {
                    900:                exthpc(vl,vx,u,i,&uh[i]); exthpc(vl,vx,v,i,&vh[i]);
                    901:        }
                    902:
                    903:        /* register degrees in each variables */
                    904:        for ( nv = 0, tvl = rvl; tvl; tvl = NEXT(tvl), nv++ );
                    905:        md = (int *)ALLOCA(nv*sizeof(int));
                    906:        for ( i = 0, tvl = rvl; i < nv; tvl = NEXT(tvl), i++ )
1.10      noro      907:                md[i] = getdeg(tvl->v,f);
1.7       noro      908:
                    909:        /* XXX for removing content of factor wrt vx */
                    910:        NEWVL(onevl); onevl->v = vx; NEXT(onevl) = 0;
                    911:
                    912:        for ( j = 1; j <= dbd; j++ ) {
                    913:                for ( i = 0, tvl = rvl; i < nv; tvl = NEXT(tvl), i++ )
                    914:                        if ( getdeg(tvl->v,u)+getdeg(tvl->v,v) > md[i] ) {
                    915:                                *up = 0;
                    916:                                return;
                    917:                        }
                    918:                for ( i = 0, t = 0; i <= j; i++ ) {
                    919:                        mulp(vl,uh[i],vh[j-i],&s); addp(vl,s,t,&w); t = w;
                    920:                }
1.10      noro      921:
1.7       noro      922:                /* s = degree j part of (f-uv) */
                    923:                exthpc(vl,vx,ff,j,&fj); subp(vl,fj,t,&s);
                    924:                for ( i = 0, wu = 0, wv = 0; i <= n; i++ ) {
1.10      noro      925:                        if ( !s )
1.7       noro      926:                                si = 0;
                    927:                        else if ( VR(s) == vx )
                    928:                                coefp(s,i,&si);
                    929:                        else if ( i == 0 )
                    930:                                si = s;
                    931:                        else
                    932:                                si = 0;
                    933:                        if ( si ) {
1.10      noro      934:                                mulp(vl,si,cv[i],&m); addp(vl,wu,m,&t); wu = t;
                    935:                                mulp(vl,si,cu[i],&m); addp(vl,wv,m,&t); wv = t;
1.7       noro      936:                        }
                    937:                }
                    938:                if ( !wu ) {
1.10      noro      939:                        gfsn_poly_to_poly(vl,u,vy,&t);
1.12      noro      940:                        shift_sf(vl,rvl,t,mev,-1,&s);
                    941:                        if ( divtp(vl,fin,s,&q) ) {
                    942:                                cont_pp_mv_sf(vl,onevl,s,&cont,up);
1.7       noro      943:                                return;
                    944:                        }
                    945:                }
                    946:                if ( !wv ) {
1.10      noro      947:                        gfsn_poly_to_poly(vl,v,vy,&t);
1.12      noro      948:                        shift_sf(vl,rvl,t,mev,-1,&s);
                    949:                        if ( divtp(vl,fin,s,&q) ) {
1.7       noro      950:                                cont_pp_mv_sf(vl,onevl,q,&cont,up);
                    951:                                return;
                    952:                        }
                    953:                }
                    954:                addp(vl,u,wu,&t); u = t;
                    955:                addp(vl,uh[j],wu,&t); uh[j] = t;
                    956:                addp(vl,v,wv,&t); v = t;
                    957:                addp(vl,vh[j],wv,&t); vh[j] = t;
                    958:        }
                    959: }
                    960:
1.12      noro      961: void adjust_coef_sf(VL vl,VL rvl,P lcu,P u0,int *mev,P *r)
1.7       noro      962: {
                    963:        P lcu0,cu;
                    964:        DCP dc0,dcu,dc;
                    965:
1.12      noro      966:        substvp_sf(vl,rvl,lcu,mev,&lcu0);
1.7       noro      967:        divsp(vl,lcu0,LC(u0),&cu);
                    968:        for ( dc0 = 0, dcu = DC(u0); dcu; dcu = NEXT(dcu) ) {
                    969:                if ( !dc0 ) {
                    970:                        NEXTDC(dc0,dc);
                    971:                        COEF(dc) = lcu;
                    972:                } else {
                    973:                        NEXTDC(dc0,dc);
                    974:                        mulp(vl,cu,COEF(dcu),&COEF(dc));
                    975:                }
                    976:                DEG(dc) = DEG(dcu);
                    977:        }
                    978:        NEXT(dc) = 0;
                    979:        MKP(VR(u0),dc0,*r);
                    980: }
                    981:
1.8       noro      982: void extended_gcd_modyk(P u0,P v0,V x,V y,int dy,P *cu,P *cv)
1.6       noro      983: {
1.8       noro      984:        BM g,h,a,b;
                    985:
                    986:        gfsn_univariate_to_sfbm(u0,dy,&g);
                    987:        gfsn_univariate_to_sfbm(v0,dy,&h);
                    988:        sfexgcd_by_hensel(g,h,dy,&a,&b);
                    989:        sfbm_to_gfsn_univariate(a,x,y,cu);
                    990:        sfbm_to_gfsn_univariate(b,x,y,cv);
                    991: }
                    992:
                    993: /* (F[y])[x] -> F[x][y] */
                    994:
                    995: void gfsn_univariate_to_sfbm(P f,int dy,BM *r)
                    996: {
                    997:        int dx,d,i;
                    998:        BM b;
                    999:        UM cy;
                   1000:        DCP dc;
                   1001:
                   1002:        dx = getdeg(VR(f),f);
                   1003:        b = BMALLOC(dx,dy);
                   1004:        DEG(b) = dy;
                   1005:        for ( dc = DC(f); dc; dc = NEXT(dc) ) {
                   1006:                /* d : degree in x, cy : poly in y */
                   1007:                d = QTOS(DEG(dc));
                   1008:                cy = BDY((GFSN)COEF(dc));
                   1009:                for ( i = DEG(cy); i >= 0; i-- )
                   1010:                        COEF(COEF(b)[i])[d] = COEF(cy)[i];
                   1011:        }
1.9       noro     1012:        for ( i = 0; i <= dy; i++ )
                   1013:                degum(COEF(b)[i],dx);
1.8       noro     1014:        *r = b;
                   1015: }
                   1016:
                   1017: void sfbm_to_gfsn_univariate(BM f,V x,V y,P *r)
                   1018: {
                   1019:        P g;
                   1020:        VL vl;
                   1021:
                   1022:        sfbmtop(f,x,y,&g);
                   1023:        NEWVL(vl); vl->v = x;
                   1024:        NEWVL(NEXT(vl)); NEXT(vl)->v = y;
                   1025:        NEXT(NEXT(vl)) = 0;
                   1026:        poly_to_gfsn_poly(vl,g,y,r);
1.6       noro     1027: }
                   1028:
1.7       noro     1029: void poly_to_gfsn_poly(VL vl,P f,V v,P *r)
1.6       noro     1030: {
1.8       noro     1031:        VL tvl,nvl0,nvl;
                   1032:        P g;
                   1033:
                   1034:        /* (x,y,...,v,...) -> (x,y,...,v) */
                   1035:        for ( nvl0 = 0, tvl = vl; tvl; tvl = NEXT(tvl) ) {
                   1036:                if ( tvl->v != v ) {
                   1037:                        NEXTVL(nvl0,nvl);
                   1038:                        nvl->v = tvl->v;
                   1039:                }
                   1040:        }
                   1041:        NEXTVL(nvl0,nvl);
                   1042:        nvl->v = v;
                   1043:        NEXT(nvl) = 0;
                   1044:        reorderp(nvl0,vl,f,&g);
                   1045:        poly_to_gfsn_poly_main(g,v,r);
                   1046: }
                   1047:
                   1048: void poly_to_gfsn_poly_main(P f,V v,P *r)
                   1049: {
                   1050:        int d;
                   1051:        UM u;
                   1052:        GFSN g;
                   1053:        DCP dc,dct,dc0;
                   1054:
1.9       noro     1055:        if ( !f )
1.8       noro     1056:                *r = f;
1.9       noro     1057:        else if ( NUM(f) || VR(f) == v ) {
1.8       noro     1058:                d = getdeg(v,f);
                   1059:                u = UMALLOC(d);
                   1060:                ptosfum(f,u);
                   1061:                MKGFSN(u,g);
                   1062:                *r = (P)g;
                   1063:        } else {
                   1064:                for ( dc0 = 0, dct = DC(f); dct; dct = NEXT(dct) ) {
                   1065:                        NEXTDC(dc0,dc);
                   1066:                        DEG(dc) = DEG(dct);
                   1067:                        poly_to_gfsn_poly_main(COEF(dct),v,&COEF(dc));
                   1068:                }
                   1069:                NEXT(dc) = 0;
                   1070:                MKP(VR(f),dc0,*r);
                   1071:        }
1.6       noro     1072: }
                   1073:
1.7       noro     1074: void gfsn_poly_to_poly(VL vl,P f,V v,P *r)
1.6       noro     1075: {
1.8       noro     1076:        VL tvl,nvl0,nvl;
                   1077:        P g;
                   1078:
                   1079:        gfsn_poly_to_poly_main(f,v,&g);
                   1080:        /* (x,y,...,v,...) -> (x,y,...,v) */
                   1081:        for ( nvl0 = 0, tvl = vl; tvl; tvl = NEXT(tvl) ) {
                   1082:                if ( tvl->v != v ) {
                   1083:                        NEXTVL(nvl0,nvl);
                   1084:                        nvl->v = tvl->v;
                   1085:                }
                   1086:        }
                   1087:        NEXTVL(nvl0,nvl);
                   1088:        nvl->v = v;
                   1089:        NEXT(nvl) = 0;
                   1090:        reorderp(vl,nvl0,g,r);
                   1091: }
                   1092:
                   1093: void gfsn_poly_to_poly_main(P f,V v,P *r)
                   1094: {
                   1095:        DCP dc,dc0,dct;
                   1096:
                   1097:        if ( !f )
                   1098:                *r = f;
                   1099:        else if ( NUM(f) ) {
                   1100:                if ( NID((Num)f) == N_GFSN )
                   1101:                        sfumtop(v,BDY((GFSN)f),r);
                   1102:                else
                   1103:                        *r = f;
                   1104:        } else {
                   1105:                for ( dc0 = 0, dct = DC(f); dct; dct = NEXT(dct) ) {
                   1106:                        NEXTDC(dc0,dc);
                   1107:                        DEG(dc) = DEG(dct);
                   1108:                        gfsn_poly_to_poly_main(COEF(dct),v,&COEF(dc));
                   1109:                }
                   1110:                NEXT(dc) = 0;
                   1111:                MKP(VR(f),dc0,*r);
                   1112:        }
1.1       noro     1113: }
1.9       noro     1114:
                   1115: void printsfum(UM f)
                   1116: {
                   1117:        int i;
                   1118:
                   1119:        for ( i = DEG(f); i >= 0; i-- ) {
                   1120:                printf("+(");
                   1121:                printf("%d",IFTOF(COEF(f)[i]));
                   1122:                printf(")*y^%d",i);
                   1123:        }
                   1124: }
                   1125:
                   1126: void printsfbm(BM f)
                   1127: {
                   1128:        int i;
                   1129:
                   1130:        for ( i = DEG(f); i >= 0; i-- ) {
                   1131:                printf("+(");
                   1132:                printsfum(COEF(f)[i]);
                   1133:                printf(")*y^%d",i);
                   1134:        }
                   1135: }
                   1136:

FreeBSD-CVSweb <freebsd-cvsweb@FreeBSD.org>