[BACK]Return to new_pd.rr CVS log [TXT][DIR] Up to [local] / OpenXM / src / asir-contrib / testing / noro

Annotation of OpenXM/src/asir-contrib/testing/noro/new_pd.rr, Revision 1.11

1.11    ! noro        1: /* $OpenXM: OpenXM/src/asir-contrib/testing/noro/new_pd.rr,v 1.10 2014/09/19 00:48:50 ohara Exp $ */
1.1       noro        2: import("gr")$
                      3: module noro_pd$
1.8       noro        4: static GBCheck,F4,EProcs,Procs,SatHomo,GBRat,SuccSat,RepColon$
1.1       noro        5:
1.8       noro        6: localf radical_membership_sat$
1.4       noro        7: localf witness$
                      8: localf get_lc,tomonic,aa,ideal_intersection_m,redbase$
1.8       noro        9: localf para_exec,nd_gr_rat,competitive_exec,call_func,call_func_serial$
1.1       noro       10: localf call_ideal_list_intersection$
1.3       noro       11: localf call_colon,call_prime_dec$
1.7       noro       12: localf prime_dec2, prime_dec_main2$
1.1       noro       13: localf first_second$
                     14: localf third$
                     15: localf locsat,iso_comp_para,extract_qj,colon_prime_dec,extract_comp$
                     16: localf separator$
1.4       noro       17: localf member,mingen,compute_gbsyz,redcoef,recompute_trace,dtop,topnum$
1.1       noro       18: localf prepost$
                     19: localf monodec0,monodec,prod$
                     20: localf extract_qd,primary_check$
                     21: localf second$
1.8       noro       22: localf gbrat,succsat,repcolon,comp_third_tdeg,comp_tord$
1.1       noro       23: localf power$
                     24:
                     25: localf syci_dec, syc_dec$
                     26: localf syca_dec,syc0_dec$
                     27:
                     28: localf find_si0,find_si1,find_si2$
                     29: localf find_ssi0,find_ssi1,find_ssi2$
                     30:
                     31: localf init_pprocs, init_eprocs, init_procs, kill_procs$
                     32:
                     33: localf sy_dec, pseudo_dec, iso_comp, prima_dec$
                     34:
                     35: localf prime_dec, prime_dec_main, lex_predec1, zprimedec, zprimadec$
                     36: localf complete_qdecomp, partial_qdecomp, partial_qdecomp0, complete_decomp$
                     37: localf partial_decomp, partial_decomp0, zprimacomp, zprimecomp$
                     38: localf fast_gb, incremental_gb, elim_gb, ldim, make_mod_subst$
                     39: localf rsgn, find_npos, gen_minipoly, indepset$
1.8       noro       40: localf maxindep, maxindep2, contraction, contraction_m, ideal_list_intersection, ideal_intersection$
1.1       noro       41: localf radical_membership, modular_radical_membership$
                     42: localf radical_membership_rep, ideal_product, saturation$
1.11    ! noro       43: localf sat, satind, sat_ind, sat_ind_var, colon, isat$
1.1       noro       44: localf ideal_colon, ideal_sat, ideal_inclusion, qd_simp_comp, qd_remove_redundant_comp$
1.7       noro       45: localf pd_simp_comp, remove_identical_comp$
1.1       noro       46: localf pd_remove_redundant_comp, ppart, sq, gen_fctr, gen_nf, gen_gb_comp$
                     47: localf gen_mptop, lcfactor, compute_deg0, compute_deg, member$
                     48: localf elimination, setintersection, setminus, sep_list$
                     49: localf first, comp_tdeg, comp_tdeg_first, tdeg, comp_by_ord, comp_by_second$
                     50: localf gbcheck,f4,sathomo,qd_check,qdb_check$
                     51:
                     52: SatHomo=0$
                     53: GBCheck=1$
                     54: GBRat=0$
1.8       noro       55: SuccSat=0$
                     56: RepColon=0$
1.1       noro       57:
                     58: #define MAX(a,b) ((a)>(b)?(a):(b))
                     59: #define ACCUM_TIME(C,R) {T1 = time(); C += (T1[0]-T0[0])+(T1[1]-T0[1]); R += (T1[3]-T0[3]); }
                     60:
                     61: def gbrat(A)
                     62: {
1.8       noro       63:        GBRat = A;
                     64: }
                     65:
                     66: def succsat(A)
                     67: {
                     68:        SuccSat = A;
                     69: }
                     70:
                     71: def repcolon(A)
                     72: {
                     73:        RepColon = A;
1.1       noro       74: }
                     75:
                     76: def gbcheck(A)
                     77: {
                     78:        if ( A ) GBCheck = 1;
                     79:        else GBCheck = -1;
                     80: }
                     81:
                     82: def f4(A)
                     83: {
                     84:        if ( A ) F4 = 1;
                     85:        else F4 = 0;
                     86: }
                     87:
                     88: def sathomo(A)
                     89: {
                     90:        if ( A ) SatHomo = 1;
                     91:        else SatHomo = 0;
                     92: }
                     93:
                     94: def init_eprocs()
                     95: {
                     96:        if ( type(NoX=getopt(nox)) == -1 ) NoX = 0;
                     97:        if ( !EProcs ) {
                     98:                if ( NoX ) {
                     99:                        P0 = ox_launch_nox();
                    100:                        P1 = ox_launch_nox();
                    101:                } else {
                    102:                        P0 = ox_launch();
                    103:                        P1 = ox_launch();
                    104:                }
                    105:                EProcs = [P0,P1];
                    106:        }
                    107: }
                    108:
                    109: def init_pprocs(N)
                    110: {
                    111:        if ( type(NoX=getopt(nox)) == -1 ) NoX = 0;
                    112:        for ( R = [], I = 0; I < N; I++ ) {
                    113:                P = NoX ? ox_launch_nox() : ox_launch();
                    114:                R = cons(P,R);
                    115:        }
                    116:        return reverse(R);
                    117: }
                    118:
                    119: def init_procs()
                    120: {
                    121:        if ( type(NoX=getopt(nox)) == -1 ) NoX = 0;
                    122:        if ( !Procs ) {
                    123:                if ( NoX ) {
                    124:                        P0 = ox_launch_nox();
                    125:                        P1 = ox_launch_nox();
                    126:                } else {
                    127:                        P0 = ox_launch();
                    128:                        P1 = ox_launch();
                    129:                }
                    130:                Procs = [P0,P1];
                    131:        }
                    132: }
                    133:
                    134: def kill_procs()
                    135: {
                    136:        if ( Procs ) {
                    137:                ox_shutdown(Procs[0]);
                    138:                ox_shutdown(Procs[1]);
                    139:                Procs = 0;
                    140:        }
                    141:        if ( EProcs ) {
                    142:                ox_shutdown(EProcs[0]);
                    143:                ox_shutdown(EProcs[1]);
                    144:                EProcs = 0;
                    145:        }
                    146: }
                    147:
                    148: def qd_check(B,V,QD)
                    149: {
                    150:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                    151:        G = nd_gr(B,V,Mod,0);
                    152:        Iso = ideal_list_intersection(map(first,QD[0]),V,0|mod=Mod);
                    153:        Emb = ideal_list_intersection(map(first,QD[1]),V,0|mod=Mod);
                    154:        GG = ideal_intersection(Iso,Emb,V,0|mod=Mod);
                    155:        return gen_gb_comp(G,GG,Mod);
                    156: }
                    157:
                    158: def primary_check(B,V)
                    159: {
                    160:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                    161:        G = nd_gr(B,V,Mod,0);
                    162:        PL = prime_dec(G,V|indep=1,mod=Mod);
                    163:        if ( length(PL) > 1 ) return 0;
                    164:        P = PL[0][0]; Y = PL[0][1];
                    165:        Z = setminus(V,Y);
                    166:        H = elim_gb(G,Z,Y,Mod,[[0,length(Z)],[0,length(Y)]]);
                    167:        H = contraction(H,Z|mod=Mod);
                    168:        H = nd_gr(H,V,Mod,0);
                    169:        if ( gen_gb_comp(G,H,Mod) ) return nd_gr(P,V,Mod,0);
                    170:        else return 0;
                    171: }
                    172:
                    173: def qdb_check(B,V,QD)
                    174: {
                    175:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                    176:        G = nd_gr(B,V,Mod,0);
                    177:        N = length(QD);
                    178:        for ( I = 0, Q = [1]; I < N; I++ )
                    179:                for ( J = 0, QL = map(first,QD[I]), L = length(QL);
                    180:                        J < L; J++ )
1.8       noro      181:                        Q = ideal_intersection_m(Q,QL[J],V,0|mod=Mod);
                    182:        Q = nd_gr(Q,V,0,0);
1.1       noro      183:        if ( !gen_gb_comp(G,Q,Mod) )
                    184:                return 0;
                    185:        for ( I = 0; I < N; I++ ) {
                    186:                T = QD[I];
                    187:                M = length(T);
                    188:                for ( J = 0; J < M; J++ ) {
                    189:                        P = primary_check(T[J][0],V|mod=Mod);
                    190:                        if ( !P ) return 0;
                    191:                        PP = nd_gr(T[J][1],V,Mod,0);
                    192:                        if ( !gen_gb_comp(P,PP,Mod) ) return 0;
                    193:                }
                    194:        }
                    195:        return 1;
                    196: }
                    197:
                    198: def extract_qd(QD,V,Ind)
                    199: {
                    200:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                    201:        N = length(Ind);
                    202:        for ( I = 0, Q = [1]; I < N; I++ )
                    203:                for ( J = 0, QL = map(first,QD[Ind[I]]), L = length(QL);
                    204:                        J < L; J++ )
                    205:                        Q = ideal_intersection(Q,QL[J],V,0|mod=Mod);
                    206:        return Q;
                    207: }
                    208:
                    209: /* SYC primary decomositions */
                    210:
                    211: def syc_dec(B,V)
                    212: {
                    213:        if ( type(SI=getopt(si)) == -1 ) SI = 2;
1.10      ohara     214:        SIFList=[noro_pd.find_ssi0, noro_pd.find_ssi1, noro_pd.find_ssi2];
1.1       noro      215:        if ( SI<0 || SI>2 ) error("sycb_dec : si should be 0,1,2");
                    216:        SIF = SIFList[SI];
                    217:
                    218:        if ( type(MaxLevel=getopt(level)) == -1 ) MaxLevel = -1;
                    219:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                    220:        if ( type(Lexdec=getopt(lexdec)) == -1 ) Lexdec = 0;
                    221:        if ( type(Time=getopt(time)) == -1 ) Time = 0;
                    222:        if ( type(Iso=getopt(iso)) == -1 ) Iso = 0;
                    223:        if ( type(Colon=getopt(colon)) == -1 ) Colon = 1;
                    224:        Ord = 0;
                    225:        Tall = time();
                    226:        C = Gt = G = fast_gb(B,V,Mod,Ord|trace=1);
                    227:        Q = []; IntQ = [1]; First = 1;
                    228:        Tpd = Tiso = Tsep = 0;
                    229:        RTpd = RTiso = RTsep = 0;
                    230:        for ( Level = 0; MaxLevel < 0 || Level <= MaxLevel; Level++ ) {
                    231:                if ( type(Gt[0])==1 ) break;
                    232: T3 = T2 = T1 = T0 = time();
                    233:                if ( First ) {
                    234:                        PtR = prime_dec(C,V|indep=1,lexdec=Lexdec,mod=Mod,radical=1);
                    235:                        Pt = PtR[0]; IntPt = PtR[1];
                    236:                        if ( gen_gb_comp(Gt,IntPt,Mod) ) {
                    237:                                /* Gt is radical and Gt = cap Pt */
                    238:                                for ( T = Pt, Qt = []; T != []; T = cdr(T) )
                    239:                                        Qt = cons([car(T)[0],car(T)[0]],Qt);
                    240:                                return append(Q,[Qt]);
                    241:                        }
                    242:                }
                    243: T1 = time(); Tpd += (T1[0]-T0[0])+(T1[1]-T0[1]); RTpd += (T1[3]-T0[3]);
                    244:                Qt = iso_comp(Gt,Pt,V,Ord|mod=Mod,first=First,iso=Iso);
                    245:                Q = append(Q,[Qt]);
                    246:                for ( T = Qt; T != []; T = cdr(T) )
                    247:                        IntQ = ideal_intersection(IntQ,car(T)[0],V,Ord
                    248:                                |mod=Mod,
                    249:                                 gbblock=[[0,length(IntQ)],[length(IntQ),length(car(T)[0])]]);
                    250:                if ( First ) { IntP = IntPt; First = 0; }
                    251:                if ( gen_gb_comp(IntQ,G,Mod) ) break;
                    252:
                    253:                M = mingen(IntQ,V);
                    254:                for ( Pt = [], C = [1], T = M; T != []; T = cdr(T) ) {
                    255:                        Ci = colon(G,car(T),V|isgb=1);
                    256:                        if ( type(Ci[0]) != 1 ) {
                    257:                                Pi = prime_dec(Ci,V|indep=1,lexdec=Lexdec,radical=1,mod=Mod);
                    258:                                C = ideal_intersection(C,Pi[1],V,Ord);
                    259:                                Pt = append(Pt,Pi[0]);
                    260:                        }
                    261:                }
                    262:                Pt = pd_simp_comp(Pt,V|first=1,mod=Mod);
                    263:                if ( Colon ) C = ideal_colon(G,IntQ,V|mod=Mod);
                    264: T2 = time(); Tiso += (T2[0]-T1[0])+(T2[1]-T1[1]); RTiso += (T2[3]-T1[3]);
                    265:                Ok = (*SIF)(C,G,IntQ,IntP,V,Ord|mod=Mod);
                    266:                Gt = append(Ok,G);
                    267: T3 = time(); Tsep += (T3[0]-T2[0])+(T3[1]-T2[1]); RTsep += (T3[3]-T2[3]);
                    268:        }
                    269:        T4 = time(); RTall = (T4[3]-Tall[3]); Tall = (T4[0]-Tall[0])+(T3[1]-Tall[1]);
                    270:        if ( Time ) {
                    271:                print(["cpu","total",Tall,"pd",Tpd,"iso",Tiso,"sep",Tsep]);
                    272:                print(["elapsed","total",RTall,"pd",RTpd,"iso",RTiso,"sep",RTsep]);
                    273:        }
                    274:        return Q;
                    275: }
                    276:
                    277: static Tint2, RTint2$
                    278:
                    279: def syci_dec(B,V)
                    280: {
                    281:        if ( type(SI=getopt(si)) == -1 ) SI = 1;
                    282:        if ( SI<0 || SI>2 ) error("sycb_assdec : si should be 0,1,2");
                    283:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                    284:        if ( type(Lexdec=getopt(lexdec)) == -1 ) Lexdec = 0;
                    285:        if ( type(Time=getopt(time)) == -1 ) Time = 0;
                    286:        if ( type(Iso=getopt(iso)) == -1 ) Iso = 0;
                    287:        if ( type(Ass=getopt(ass)) == -1 ) Ass = 0;
                    288:        if ( type(Colon=getopt(colon)) == -1 ) Colon = 0;
                    289:        if ( type(Para=getopt(para)) == -1 ) Para = 0;
1.8       noro      290:        if ( type(Trace=getopt(trace)) == -1 ) Trace = 0;
1.1       noro      291:        Ord = 0;
                    292:        Tiso = Tint = Tpd = Text = Tint2 = 0;
                    293:        RTiso = RTint = RTpd = RText = RTint2 = 0;
                    294: T00 = time();
1.8       noro      295:        G = fast_gb(B,V,Mod,Ord|trace=Trace);
1.1       noro      296:        IntQ = [1]; QL = RL = []; First = 1;
                    297:        for ( Level = 0; ; Level++ ) {
                    298: T0 = time();
1.8       noro      299:                if ( !Level ) {
1.1       noro      300:                        PtR = prime_dec(G,V|indep=1,lexdec=Lexdec,mod=Mod,radical=1);
1.8       noro      301: ACCUM_TIME(Tfpd,RTfpd)
1.1       noro      302:                        Pt = PtR[0]; IntPt = PtR[1]; Rad = IntPt;
1.2       noro      303:                        if ( gen_gb_comp(G,Rad,Mod) ) {
                    304:                                /* Gt is radical and Gt = cap Pt */
                    305:                                for ( T = Pt, Qt = []; T != []; T = cdr(T) )
                    306:                                        Qt = cons([car(T)[0],car(T)[0],car(T)[1]],Qt);
                    307:                                return [reverse(Qt)];
                    308:                        }
1.1       noro      309:                } else
                    310:                        Pt = colon_prime_dec(G,IntQ,V|lexdec=Lexdec,mod=Mod,para=Para);
                    311: ACCUM_TIME(Tpd,RTpd)
                    312: T0 = time();
                    313:                Rt = iso_comp(G,Pt,V,Ord|mod=Mod,iso=Iso,para=Para,intq=IntQ);
                    314: ACCUM_TIME(Tiso,RTiso)
1.8       noro      315:                if ( !Level ) {
                    316:                        if ( Iso == 3 ) {
                    317:                                NI = length(Rt);
                    318:                                Q = IntQ;
                    319: T0 = time();
                    320:                                if ( Para ) {
                    321:                                        for ( J = 0, Task = []; J < NI; J++ ) {
                    322:                                                T = ["noro_pd.extract_qj",Q,V,Rt[J],Rad,Mod,SI,Colon,-1];
                    323:                                                Task = cons(T,Task);
                    324:                                        }
                    325:                                        Task = reverse(Task);
                    326:                                        print("comps:",2); print(length(Task),2);
                    327:                                        Rt = para_exec(Para,Task);
                    328:                                } else {
                    329:                                        for ( J = 0, T = []; J < NI; J++ ) {
                    330:                                                TJ = extract_qj(Q,V,Rt[J],Rad,Mod,SI,Colon,-1);
                    331:                                                T = cons(TJ,T);
                    332:                                        }
                    333:                                        Rt = reverse(T);
                    334:                                }
                    335: ACCUM_TIME(Text,RText)
                    336:                        }
                    337:                        print("");
                    338: T0 = time();
                    339:                        Int = Rad;
                    340:                        for ( T = Rt; T != []; T = cdr(T) )
                    341:                                if ( !gb_comp(car(T)[0],car(T)[1]) )
                    342:                                        Int = ideal_intersection_m(Int,car(T)[0],V,Ord|mod=Mod);
                    343:                        IntQ = nd_gr(Int,V,Mod,Ord);
                    344: ACCUM_TIME(Tint,RTint)
                    345:                        RL = append(RL,[Rt]);
                    346:                } else if ( Iso != 3 ) {
1.1       noro      347: T0 = time();
1.8       noro      348:                        IntQ = ideal_list_intersection(map(first,Rt),V,Ord|mod=Mod,isgb=1);
1.4       noro      349:                        RL = append(RL,[Rt]);
1.8       noro      350: ACCUM_TIME(Tint,RTint)
1.4       noro      351:                } else {
                    352:                        NI = length(Rt);
                    353:                        Q = IntQ;
1.8       noro      354:                        if ( Para ) {
                    355:                                for ( J = 0, Task = []; J < NI; J++ ) {
                    356:                                        T = ["noro_pd.extract_qj",Q,V,Rt[J],Rad,Mod,SI,Colon,-1];
                    357:                                        Task = cons(T,Task);
                    358:                                }
                    359:                                Task = reverse(Task);
                    360:                                print("comps:",2); print(length(Task),2);
                    361: T0 = time();
                    362:                                R = para_exec(Para,Task);
                    363: ACCUM_TIME(Text,RText)
                    364:                                print("");
                    365: T0 = time();
                    366:                                IntQ = ideal_list_intersection(cons(IntQ,map(first,R)),V,Ord|mod=Mod);
                    367: ACCUM_TIME(Tint,RTint)
                    368:                                RL = append(RL,[R]);
                    369:                        } else {
                    370:                                for ( J = 0, T = []; J < NI; J++ ) {
                    371: T0 = time();
                    372:                                        TJ = extract_qj(Q,V,Rt[J],Rad,Mod,SI,Colon,-1);
                    373: ACCUM_TIME(Text,RText)
                    374:                                        T = cons(TJ,T);
                    375: T0 = time();
                    376:                                        IntQ = ideal_intersection_m(IntQ,TJ[0],V,Ord|mod=Mod);
                    377: ACCUM_TIME(Tint,RTint)
                    378:                                }
                    379:                                print("");
                    380: T0 = time();
                    381:                                IntQ = nd_gr(IntQ,V,Mod,Ord);
                    382: ACCUM_TIME(Tint,RTint)
                    383:                                T = reverse(T); RL = append(RL,[T]);
1.4       noro      384:                        }
                    385:                }
1.1       noro      386:                QL = append(QL,[IntQ]);
                    387:                if ( gen_gb_comp(IntQ,G,Mod) ) break;
                    388:        }
                    389: T0 = time();
1.4       noro      390:        if ( Iso != 3 && !Ass )
1.1       noro      391:                RL = extract_comp(QL,RL,V,Rad|mod=Mod,para=Para,si=SI,colon=Colon,ass=Ass);
                    392: ACCUM_TIME(Text,RText)
                    393:        if ( Time ) {
                    394:                T1 = time();
                    395:                Tall = T1[0]-T00[0]+T1[1]-T00[1]; RTall += T1[3]-T00[3];
                    396:                Tass = Tall-Text; RTass = RTall-RText;
1.8       noro      397:                print(["total",Tall,"ass",Tass,"pd",Tpd,"(fpd)",Tfpd,"iso",Tiso,"int",Tint,"ext",Text]);
                    398:                print(["elapsed",RTall,"ass",RTass,"pd",RTpd,"(fpd)",RTfpd,"iso",RTiso,"int",RTint,"ext",RText]);
1.1       noro      399:        }
                    400:        return  RL;
                    401: }
                    402:
                    403: def extract_comp(QL,RL,V,Rad) {
                    404:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                    405:        if ( type(Para=getopt(para)) == -1 ) Para = 0;
                    406:        if ( type(Colon=getopt(colon)) == -1 ) Colon = 0;
                    407:        if ( type(SI=getopt(si)) == -1 ) SI = 1;
                    408:        if ( type(Ass=getopt(ass)) == -1 ) Ass = 0;
                    409:
                    410:        L = length(QL);
                    411:        if ( Para ) {
                    412:                for ( Task = [], I = 1; I < L; I++ ) {
                    413:                        QI = QL[I-1]; RI = RL[I]; NI = length(RI);
                    414:                        for ( J = 0, TI = []; J < NI; J++ ) {
                    415:                                T = ["noro_pd.extract_qj",QI,V,RI[J],Rad,Mod,SI,Colon,I];
                    416:                                Task = cons(T,Task);
                    417:                        }
                    418:                }
1.8       noro      419:                Task = reverse(Task);
1.1       noro      420:                print("comps:",2); print(length(Task),2); print("");
                    421:                R = para_exec(Para,Task);
                    422:                S = vector(L);
                    423:                for ( I = 1; I < L; I++ ) S[I] = [];
                    424:                S[0] = RL[0];
                    425:                for ( T = R; T != []; T = cdr(T) ) {
                    426:                        U = car(T); Level = U[0]; Body = U[1];
                    427:                        S[Level] = cons(Body,S[Level]);
                    428:                }
                    429:                return vtol(S);
                    430:        } else {
                    431:                TL = [RL[0]];
                    432:                for ( I = 1; I < L; I++ ) {
                    433:                        print("level:",2); print(I,2);
                    434:                        print(" comps:",2); print(length(RL[I]),2); print("");
                    435:                        QI = QL[I-1]; RI = RL[I]; NI = length(RI);
                    436:                        for ( J = 0, TI = []; J < NI; J++ ) {
                    437:                                TIJ = extract_qj(QI,V,RI[J],Rad,Mod,SI,Colon,-1);
                    438:                                TI = cons(TIJ,TI);
                    439:                        }
                    440:                        TI = reverse(TI); TL = cons(TI,TL);
                    441:                }
                    442:                TL = reverse(TL);
                    443:        }
                    444:        return TL;
                    445: }
                    446:
                    447: def colon_prime_dec(G,IntQ,V) {
                    448:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                    449:        if ( type(Lexdec=getopt(lexdec)) == -1 ) Lexdec = 0;
                    450:        if ( type(Para=getopt(para)) == -1 ) Para = 0;
1.3       noro      451:        if ( !Para ) {
                    452:                print("colon_pd:",2); print(length(IntQ),2);
                    453:        }
1.1       noro      454:        if ( !Mod ) M = mingen(IntQ,V);
                    455:        else M = IntQ;
                    456:        if ( Para ) {
                    457:                L = length(M);
1.3       noro      458:                for ( Task = [], J = 0; J < L; J++ )
1.1       noro      459:                        if ( gen_nf(M[J],G,V,Ord,Mod) ) {
1.3       noro      460:                                T = ["noro_pd.call_colon",G,M[J],V,Mod,1];
1.1       noro      461:                                Task = cons(T,Task);
                    462:                        }
                    463:                Task = reverse(Task);
                    464:                R = para_exec(Para,Task);
1.3       noro      465:                R = pd_simp_comp(R,V|mod=Mod); L = length(R);
                    466:
                    467:                for ( Task = [], J = 0; J < L; J++ ) {
                    468:                        T = ["noro_pd.call_prime_dec",R[J],V,1,Lexdec,Mod];
                    469:                        Task = cons(T,Task);
                    470:                }
                    471:                Task = reverse(Task);
                    472:                R = para_exec(Para,Task);
                    473:
1.1       noro      474:                for ( Pt = [], T = R; T != []; T = cdr(T) ) Pt = append(Pt,car(T));
                    475:        } else {
1.3       noro      476:                for ( R = [], T = M; T != []; T = cdr(T) ) {
                    477:                        Ci = colon(G,car(T),V|isgb=1,mod=Mod);
                    478:                        R = cons(Ci,R);
                    479:                }
                    480:                print("->",2); print(length(M),2);
                    481:                R = pd_simp_comp(R,V|mod=Mod);
                    482:                print("->",2); print(length(R));
1.8       noro      483: #if 1
1.3       noro      484:                for ( Pt = [], T = R; T != []; T = cdr(T) ) {
                    485:                        Pi = prime_dec(car(T),V|indep=1,lexdec=Lexdec,mod=Mod);
1.1       noro      486:                        Pt = append(Pt,Pi);
                    487:                }
1.8       noro      488: #else
                    489:                J = ideal_list_intersection(R,V,0|mod=Mod);
                    490:                Pt = prime_dec(J,V|indep=1,lexdec=Lexdec,mod=Mod);
                    491: #endif
1.1       noro      492:        }
1.8       noro      493: #if 1
1.1       noro      494:        Pt = pd_simp_comp(Pt,V|first=1,mod=Mod);
1.8       noro      495: #endif
1.1       noro      496:        return Pt;
                    497: }
                    498:
1.3       noro      499: def call_colon(G,F,V,Mod,IsGB)
                    500: {
                    501:        return colon(G,F,V|isgb=1,mod=Mod);
                    502: }
                    503:
                    504: def call_prime_dec(G,V,Indep,Lexdec,Mod)
1.1       noro      505: {
1.3       noro      506:        if ( type(G[0]) != 1 )
                    507:                Pi = prime_dec(G,V|indep=Indep,lexdec=Lexdec,mod=Mod);
1.1       noro      508:        else
                    509:                Pi = [];
                    510:        return Pi;
                    511: }
                    512:
                    513: def extract_qj(Q,V,QL,Rad,Mod,SI,Colon,Level)
                    514: {
1.10      ohara     515:        SIFList=[noro_pd.find_ssi0, noro_pd.find_ssi1, noro_pd.find_ssi2];
1.1       noro      516:        SIF = SIFList[SI];
                    517:        G = QL[0]; P = QL[1]; PV = QL[2];
1.4       noro      518:        if ( Q != [1] ) {
                    519:                C = Colon ? ideal_colon(G,Q,V|mod=Mod) : P;
                    520:                Ok = (*SIF)(C,G,Q,Rad,V,0|mod=Mod);
                    521:        } else
                    522:                Ok = [];
1.1       noro      523:        V0 = setminus(V,PV);
                    524:        HJ = elim_gb(append(Ok,G),V0,PV,Mod,[[0,length(V0)],[0,length(PV)]]);
                    525:        HJ = contraction(HJ,V0|mod=Mod);
                    526:        IJ = nd_gr(HJ,V,Mod,Ord);
                    527:        return Level >= 0 ? [Level,[IJ,P]] : [IJ,P];
                    528: }
                    529:
                    530: def syca_dec(B,V)
                    531: {
                    532: T00 = time();
                    533:        if ( type(SI=getopt(si)) == -1 ) SI = 2;
1.10      ohara     534:        SIFList=[noro_pd.find_si0, noro_pd.find_si1, noro_pd.find_si2];
                    535:        SIF = SIFList[SI];
1.1       noro      536:        if ( !SIF ) error("syca_dec : si should be 0,1,2");
                    537:
                    538:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                    539:        if ( type(Lexdec=getopt(lexdec)) == -1 ) Lexdec = 0;
                    540:        if ( type(NoSimp=getopt(nosimp)) == -1 ) NoSimp = 0;
                    541:        if ( type(Time=getopt(time)) == -1 ) Time = 0;
                    542:        if ( type(Iso=getopt(iso)) == -1 ) Iso = 0;
                    543:        Ord = 0;
                    544:        Gt = G0 = G = fast_gb(B,V,Mod,Ord|trace=1);
                    545:        Q0 = Q = []; IntQ0 = IntQ = [1]; First = 1;
                    546:        C = 0;
                    547:
                    548:        Tass = Tiso = Tcolon = Tsep = Tirred = 0;
                    549:        Rass = Riso = Rcolon = Rsep = Rirred = 0;
                    550:        while ( 1 ) {
                    551:                if ( type(Gt[0])==1 ) break;
                    552:                T0 = time();
                    553:                PtR = prime_dec(Gt,V|indep=1,lexdec=Lexdec,mod=Mod,radical=1);
                    554:                T1 = time(); Tass += T1[0]-T0[0]+T1[1]-T0[1]; Rass += T1[3]-T0[3];
                    555:                Pt = PtR[0]; IntPt = PtR[1];
                    556:                if ( gen_gb_comp(Gt,IntPt,Mod) ) {
                    557:                        /* Gt is radical and Gt = cap Pt */
                    558:                        for ( T = Pt, Qt = []; T != []; T = cdr(T) )
                    559:                                Qt = cons([car(T)[0],car(T)[0]],Qt);
                    560:                        if ( First )
                    561:                                return [Qt,[]];
                    562:                        else
                    563:                                Q0 = append(Qt,Q0);
                    564:                        break;
                    565:                }
                    566:                T0 = time();
                    567:                Qt = iso_comp(Gt,Pt,V,Ord|mod=Mod,isgb=1,iso=Iso);
                    568:                T1 = time(); Tiso += T1[0]-T0[0]+T1[1]-T0[1]; Riso += T1[3]-T0[3];
                    569:                IntQt = ideal_list_intersection(map(first,Qt),V,Ord|mod=Mod);
                    570:                if ( First ) {
                    571:                        IntQ0 = IntQ = IntQt; IntP = IntPt; Qi = Qt; First = 0;
                    572:                } else {
                    573:                        IntQ1 = ideal_intersection(IntQ,IntQt,V,Ord|mod=Mod);
                    574:                        if ( gen_gb_comp(IntQ,IntQ1,Mod) ) {
                    575:                                G = Gt; IntP = IntPt; Q = []; IntQ = [1]; C = 0;
                    576:                                continue;
                    577:                        } else {
                    578:                                IntQ = IntQ1;
                    579:                                IntQ1 = ideal_intersection(IntQ0,IntQt,V,Ord|mod=Mod);
                    580:                                if ( !gen_gb_comp(IntQ0,IntQ1,Mod) ) {
                    581:                                        Q = append(Qt,Q);
                    582:                                        for ( T = Qt; T != []; T = cdr(T) )
                    583:                                                if ( !ideal_inclusion(IntQ0,car(T)[0],V,Ord|mod=Mod) )
                    584:                                                        Q0 = append(Q0,[car(T)]);
                    585:                                        IntQ0 = IntQ1;
                    586:                                }
                    587:                        }
                    588:                }
                    589:                if ( gen_gb_comp(IntQt,Gt,Mod) || gen_gb_comp(IntQ,G,Mod) || gen_gb_comp(IntQ0,G0,Mod) ) break;
                    590:                T0 = time();
                    591:                C1 = ideal_colon(G,IntQ,V|mod=Mod);
                    592:                T1 = time(); Tcolon += T1[0]-T0[0]+T1[1]-T0[1]; Rcolon += T1[3]-T0[3];
                    593:                if ( C && gen_gb_comp(C,C1,Mod) ) {
                    594:                        G = Gt; IntP = IntPt; Q = []; IntQ = [1]; C = 0;
                    595:                        continue;
                    596:                } else C = C1;
                    597:                T0 = time();
                    598:                Ok = (*SIF)(C,G,IntQ,IntP,V,Ord|mod=Mod);
                    599:                G1 = append(Ok,G);
                    600:                Gt1 = incremental_gb(G1,V,Ord|mod=Mod);
                    601:                T1 = time(); Tsep += T1[0]-T0[0]+T1[1]-T0[1]; Rsep += T1[3]-T0[3];
                    602:                Gt = Gt1;
                    603:        }
                    604:        T0 = time();
                    605:        if ( !NoSimp ) Q1 = qd_remove_redundant_comp(G0,Qi,Q0,V,Ord|mod=Mod);
                    606:        else Q1 = Q0;
                    607:        if ( Time ) {
                    608:                T1 = time(); Tirred += T1[0]-T0[0]+T1[1]-T0[1]; Rirred += T1[3]-T0[3];
                    609:                Tall = T1[0]-T00[0]+T1[1]-T00[1]; Rall += T1[3]-T00[3];
                    610:                print(["total",Tall,"ass",Tass,"iso",Tiso, "colon",Tcolon,"sep",Tsep,"irred",Tirred]);
                    611:                print(["Rtotal",Rall,"Rass",Rass,"Riso",Riso, "Rcolon",Rcolon,"Rsep",Rsep,"Rirred",Rirred]);
                    612:                print(["iso",length(Qi),"emb",length(Q0),"->",length(Q1)]);
                    613:        }
                    614:        return [Qi,Q1];
                    615: }
                    616:
                    617: def syc0_dec(B,V)
                    618: {
                    619: T00 = time();
                    620:        if ( type(SI=getopt(si)) == -1 ) SI = 1;
1.10      ohara     621:        SIFList=[noro_pd.find_si0, noro_pd.find_si1, noro_pd.find_si2, noro_pd.find_ssi0, noro_pd.find_ssi1, noro_pd.find_ssi2];
                    622:        SIF = SIFList[SI];
1.1       noro      623:        if ( !SIF ) error("syc0_dec : si should be 0,1,2");
                    624:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                    625:        if ( type(Lexdec=getopt(lexdec)) == -1 ) Lexdec = 0;
                    626:        if ( type(NoSimp=getopt(nosimp)) == -1 ) NoSimp = 0;
                    627:        if ( type(Time=getopt(time)) == -1 ) Time = 0;
                    628:        Ord = 0;
                    629:        G = fast_gb(B,V,Mod,Ord);
                    630:        Q = []; IntQ = [1]; Gt = G; First = 1;
                    631:        Tass = Tiso = Tcolon = Tsep = Tirred = 0;
                    632:        Rass = Riso = Rcolon = Rsep = Rirred = 0;
                    633:        while ( 1 ) {
                    634:                if ( type(Gt[0])==1 ) break;
                    635:                T0 = time();
                    636:                PtR = prime_dec(Gt,V|indep=1,lexdec=Lexdec,mod=Mod,radical=1);
                    637:                T1 = time(); Tass += T1[0]-T0[0]+T1[1]-T0[1]; Rass += T1[3]-T0[3];
                    638:                Pt = PtR[0]; IntPt = PtR[1];
                    639:                if ( gen_gb_comp(Gt,IntPt,Mod) ) {
                    640:                        /* Gt is radical and Gt = cap Pt */
                    641:                        for ( T = Pt, Qt = []; T != []; T = cdr(T) )
                    642:                                Qt = cons([car(T)[0],car(T)[0]],Qt);
                    643:                        if ( First )
                    644:                                return [Qt,[]];
                    645:                        else
                    646:                                Q = append(Qt,Q);
                    647:                        break;
                    648:                }
                    649:
                    650:                T0 = time();
                    651:                Qt = iso_comp(Gt,Pt,V,Ord|mod=Mod,isgb=1);
                    652:                T1 = time(); Tiso += T1[0]-T0[0]+T1[1]-T0[1]; Riso += T1[3]-T0[3];
                    653:                IntQt = ideal_list_intersection(map(first,Qt),V,Ord|mod=Mod);
                    654:                if ( First ) {
                    655:                        IntQ = IntQt; Qi = Qt; First = 0;
                    656:                } else {
                    657:                        IntQ1 = ideal_intersection(IntQ,IntQt,V,Ord|mod=Mod);
                    658:                        if ( !gen_gb_comp(IntQ1,IntQ,Mod) )
                    659:                                Q = append(Qt,Q);
                    660:                }
                    661:                if ( gen_gb_comp(IntQ,G,Mod) || gen_gb_comp(IntQt,Gt,Mod) )
                    662:                        break;
                    663:                T0 = time();
                    664:                C = ideal_colon(Gt,IntQt,V|mod=Mod);
                    665:                T1 = time(); Tcolon += T1[0]-T0[0]+T1[1]-T0[1]; Rcolon += T1[3]-T0[3];
                    666:                T0 = time();
                    667:                Ok = (*SIF)(C,Gt,IntQt,IntPt,V,Ord|mod=Mod);
                    668:                G1 = append(Ok,Gt);
                    669:                Gt = incremental_gb(G1,V,Ord|mod=Mod);
                    670:                T1 = time(); Tsep += T1[0]-T0[0]+T1[1]-T0[1]; Rsep += T1[3]-T0[3];
                    671:        }
                    672:        T0 = time();
                    673:        if ( !NoSimp ) Q1 = qd_remove_redundant_comp(G,Qi,Q,V,Ord|mod=Mod);
                    674:        else Q1 = Q;
                    675:        T1 = time(); Tirred += T1[0]-T0[0]+T1[1]-T0[1]; Rirred += T1[3]-T0[3];
                    676:        Tall = T1[0]-T00[0]+T1[1]-T00[1]; Rall += T1[3]-T00[3];
                    677:        if ( Time ) {
                    678:                print(["total",Tall,"ass",Tass,"iso",Tiso, "colon",Tcolon,"sep",Tsep,"irred",Tirred]);
                    679:                print(["Rtotal",Rall,"Rass",Rass,"Riso",Riso, "Rcolon",Rcolon,"Rsep",Rsep,"Rirred",Rirred]);
                    680:                print(["iso",length(Qi),"emb",length(Q),"->",length(Q1)]);
                    681:        }
                    682:        return [Qi,Q1];
                    683: }
                    684:
                    685: def power(A,I) { return A^I; }
                    686:
                    687:
                    688: /* functions for computating a separing ideal  */
                    689: /* C=G:Q, Rad=rad(Q), return J s.t. Q cap (G+J) = G */
                    690:
                    691: def find_si0(C,G,Q,Rad,V,Ord) {
                    692:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                    693:        for ( CI = C, I = 1; ; I++ ) {
                    694:                for ( T = CI, S = []; T != []; T = cdr(T) )
                    695:                        if ( gen_nf(car(T),Q,V,Ord,Mod) ) S = cons(car(T),S);
                    696:                if ( S == [] )
                    697:                        error("find_si0 : cannot happen");
                    698:                G1 = append(S,G);
                    699:                Int = ideal_intersection(G1,Q,V,Ord|mod=Mod);
                    700:                /* check whether (Q cap (G+S)) = G */
                    701:                if ( gen_gb_comp(Int,G,Mod) ) { print([0]); return reverse(S); }
                    702:                CI = ideal_product(CI,C,V|mod=Mod);
                    703:        }
                    704: }
                    705:
                    706: def find_si1(C,G,Q,Rad,V,Ord) {
                    707:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                    708:        for ( T = C, S = []; T != []; T = cdr(T) )
                    709:                if ( gen_nf(car(T),Q,V,Ord,Mod) ) S = cons(car(T),S);
                    710:        if ( S == [] )
                    711:                error("find_si1 : cannot happen");
                    712:        G1 = append(S,G);
                    713:        Int = ideal_intersection(G1,Q,V,Ord|mod=Mod);
                    714:        /* check whether (Q cap (G+S)) = G */
                    715:        if ( gen_gb_comp(Int,G,Mod) ) { print([0]); return reverse(S); }
                    716:
1.9       ohara     717:        C = qsort(C,noro_pd.comp_tdeg);
1.1       noro      718:
                    719:        Tmp = ttttt; TV = cons(Tmp,V); Ord1 = [[0,1],[Ord,length(V)]];
                    720:        Int0 = incremental_gb(append(vtol(ltov(G)*Tmp),vtol(ltov(Q)*(1-Tmp))),
                    721:                TV,Ord1|gbblock=[[0,length(G)]],mod=Mod);
                    722:        Dp = dp_gr_print(); dp_gr_print(0);
                    723:        for ( T = C, S = []; T != []; T = cdr(T) ) {
                    724:                if ( !gen_nf(car(T),Rad,V,Ord,Mod) ) continue;
                    725:                Ui = U = car(T);
                    726:                for ( I = 1; ; I++ ) {
                    727:                        G1 = cons(Ui,G);
                    728:                        Int = ideal_intersection(G1,Q,V,Ord|mod=Mod);
                    729:                        if ( gen_gb_comp(Int,G,Mod) ) break;
                    730:                        else
                    731:                                Ui = gen_nf(Ui*U,G,V,Ord,Mod);
                    732:                }
                    733:                print([length(T),I],2);
                    734:                Int1 = incremental_gb(append(Int0,[Tmp*Ui]),TV,Ord1
                    735:                        |gbblock=[[0,length(Int0)]],mod=Mod);
                    736:                Int = elimination(Int1,V);
                    737:                if ( !gen_gb_comp(Int,G,Mod) ) {
                    738:                        break;
                    739:                } else {
                    740:                        Int0 = Int1;
                    741:                        S = cons(Ui,S);
                    742:                }
                    743:        }
                    744:        print("");
                    745:        dp_gr_print(Dp);
                    746:        return reverse(S);
                    747: }
                    748:
                    749: def find_si2(C,G,Q,Rad,V,Ord) {
                    750:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                    751:        for ( T = C, S = []; T != []; T = cdr(T) )
                    752:                if ( gen_nf(car(T),Q,V,Ord,Mod) ) S = cons(car(T),S);
                    753:        if ( S == [] )
                    754:                error("find_si2 : cannot happen");
                    755:        G1 = append(S,G);
                    756:        Int = ideal_intersection(G1,Q,V,Ord|mod=Mod);
                    757:        /* check whether (Q cap (G+S)) = G */
                    758:        if ( gen_gb_comp(Int,G,Mod) ) { print([0]); return reverse(S); }
                    759:
1.9       ohara     760:        C = qsort(C,noro_pd.comp_tdeg);
1.1       noro      761:
                    762:        Dp = dp_gr_print(); dp_gr_print(0);
                    763:        Tmp = ttttt; TV = cons(Tmp,V); Ord1 = [[0,1],[Ord,length(V)]];
                    764:        Int0 = incremental_gb(append(vtol(ltov(G)*Tmp),vtol(ltov(Q)*(1-Tmp))),
                    765:                TV,Ord1|gbblock=[[0,length(G)]],mod=Mod);
                    766:        for ( T = C, S = []; T != []; T = cdr(T) ) {
                    767:                if ( !gen_nf(car(T),Rad,V,Ord,Mod) ) continue;
                    768:                Ui = U = car(T);
                    769:                for ( I = 1; ; I++ ) {
                    770:                        Int1 = incremental_gb(append(Int0,[Tmp*Ui]),TV,Ord1
                    771:                                        |gbblock=[[0,length(Int0)]],mod=Mod);
                    772:                        Int = elimination(Int1,V);
                    773:                        if ( gen_gb_comp(Int,G,Mod) ) break;
                    774:                        else
                    775:                                Ui = gen_nf(Ui*U,G,V,Ord,Mod);
                    776:                }
                    777:                print([length(T),I],2);
                    778:                S = cons(Ui,S);
                    779:        }
1.9       ohara     780:        S = qsort(S,noro_pd.comp_tdeg);
1.1       noro      781:        print("");
                    782:        End = Len = length(S);
                    783:
                    784:        Tmp = ttttt; TV = cons(Tmp,V); Ord1 = [[0,1],[Ord,length(V)]];
                    785:        Prev = 1;
                    786:        G1 = append(G,[S[0]]);
                    787:        Int0 = incremental_gb(append(vtol(ltov(G1)*Tmp),vtol(ltov(Q)*(1-Tmp))),
                    788:                TV,Ord1|gbblock=[[0,length(G)]],mod=Mod);
                    789:        if ( End > 1 ) {
                    790:                Cur = 2;
                    791:                while ( Prev < Cur ) {
                    792:                        for ( St = [], I = Prev; I < Cur; I++ ) St = cons(Tmp*S[I],St);
                    793:                        Int1 = incremental_gb(append(Int0,St),TV,Ord1
                    794:                                |gbblock=[[0,length(Int0)]],mod=Mod);
                    795:                        Int = elimination(Int1,V);
                    796:                        if ( gen_gb_comp(Int,G,Mod) ) {
                    797:                                print([Cur],2);
                    798:                                Prev = Cur;
                    799:                                Cur = Cur+idiv(End-Cur+1,2);
                    800:                                Int0 = Int1;
                    801:                        } else {
                    802:                                End = Cur;
                    803:                                Cur = Prev + idiv(Cur-Prev,2);
                    804:                        }
                    805:                }
                    806:                for ( St = [], I = 0; I < Prev; I++ ) St = cons(S[I],St);
                    807:        } else
                    808:                St = [S[0]];
                    809:        print("");
                    810:        for ( I = Prev; I < Len; I++ ) {
                    811:                Int1 = incremental_gb(append(Int0,[Tmp*S[I]]),TV,Ord1
                    812:                        |gbblock=[[0,length(Int0)]],mod=Mod);
                    813:                Int = elimination(Int1,V);
                    814:                if ( gen_gb_comp(Int,G,Mod) ) {
                    815:                        print([I],2);
                    816:                        St = cons(S[I],St);
                    817:                        Int0 = Int1;
                    818:                }
                    819:        }
                    820:        Ok = reverse(St);
                    821:        print("");
                    822:        print([length(S),length(Ok)]);
                    823:        dp_gr_print(Dp);
                    824:        return Ok;
                    825: }
                    826:
                    827: /* functions for computing a saturated separating ideal */
                    828:
                    829: def find_ssi0(C,G,Q,Rad,V,Ord) {
                    830:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                    831:        if ( type(Reduce=getopt(red)) == -1 ) Reduce = 0;
                    832:        for ( T = C, S = []; T != []; T = cdr(T) )
                    833:                if ( gen_nf(car(T),Q,V,Ord,Mod) ) S = cons(car(T),S);
                    834:        if ( S == [] )
                    835:                error("find_ssi0 : cannot happen");
                    836:        G1 = append(S,G);
                    837:        Int = ideal_intersection(G1,Q,V,Ord|mod=Mod);
                    838:        /* check whether (Q cap (G+S)) = G */
                    839:        if ( gen_gb_comp(Int,G,Mod) ) { print([0]); return reverse(S); }
                    840:
                    841:        if ( Reduce ) {
                    842:                for ( T = C, U = []; T != []; T = cdr(T) )
                    843:                        if ( gen_nf(car(T),Rad,V,Ord,Mod) ) U = cons(car(T),U);
                    844:                U = reverse(U);
                    845:        } else
                    846:                U = C;
                    847:
                    848:        for ( I = 1; ; I++ ) {
                    849:                print([I],2);
                    850:                S = map(power,U,I);
                    851:                G1 = append(S,G);
                    852:                Int = ideal_intersection(G1,Q,V,Ord|mod=Mod);
                    853:                /* check whether (Q cap (G+S)) = G */
                    854:                if ( gen_gb_comp(Int,G,Mod) ) { print(""); return reverse(S); }
                    855:        }
                    856: }
                    857:
                    858: def find_ssi1(C,G,Q,Rad,V,Ord) {
                    859:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                    860:        if ( type(Reduce=getopt(red)) == -1 ) Reduce = 0;
                    861:        for ( T = C, S = []; T != []; T = cdr(T) )
                    862:                if ( gen_nf(car(T),Q,V,Ord,Mod) ) S = cons(car(T),S);
                    863:        if ( S == [] )
                    864:                error("find_ssi1 : cannot happen");
                    865:        G1 = append(S,G);
                    866:        Int = ideal_intersection(G1,Q,V,Ord|mod=Mod);
                    867:        /* check whether (Q cap (G+S)) = G */
                    868:        if ( gen_gb_comp(Int,G,Mod) ) { print([0]); return reverse(S); }
                    869:
                    870:        dp_ord(Ord); DC = map(dp_ptod,C,V);
1.9       ohara     871:        DC = qsort(DC,noro_pd.comp_tord); C = map(dp_dtop,DC,V);
1.1       noro      872:        print(length(C),2);
                    873:        if ( Reduce ) {
                    874:                SC = map(sq,C,Mod);
                    875:                SC = reverse(SC); C = reverse(C);
                    876:                for ( T = C, C1 = [],  R1 = append(SC,Rad); T != []; T = cdr(T) ) {
                    877:                        R0 = car(R1); R1 = cdr(R1);
                    878:                        if ( !gen_nf(car(T),Rad,V,Ord,Mod) ) continue;
                    879:                        if ( radical_membership(R0,R1,V|mod=Mod)  ) {
                    880:                                C1 = cons(car(T),C1);
                    881:                                R1 = append(R1,[R0]);
                    882:                        }
                    883:                }
                    884:                print("->",0); print(length(C1),2);
                    885:                C = C1;
                    886:        }
                    887:        print(" ",2);
                    888:
                    889:        Tmp = ttttt; TV = cons(Tmp,V); Ord1 = [[0,1],[Ord,length(V)]];
                    890:        Int0 = incremental_gb(append(vtol(ltov(G)*Tmp),vtol(ltov(Q)*(1-Tmp))),
                    891:                TV,Ord1|gbblock=[[0,length(G)]],mod=Mod);
                    892:        Dp = dp_gr_print(); dp_gr_print(0);
                    893:        for ( J = 0, T = C, S = [], GS = G; T != []; T = cdr(T), J++ ) {
                    894:                if ( !gen_nf(car(T),Rad,V,Ord,Mod) ) continue;
                    895:                Ui = U = car(T);
                    896:                for ( I = 1; ; I++ ) {
                    897:                        Int1 = nd_gr(append(Int0,[Tmp*Ui]),TV,Mod,Ord1
                    898:                                |gbblock=[[0,length(Int0)]],newelim=1);
                    899:                        if ( Int1 ) {
                    900:                                Int = elimination(Int1,V);
                    901:                                if ( gen_gb_comp(Int,G,Mod) ) break;
                    902:                        }
                    903:                        print("x",2);
                    904:                        Ui = gen_nf(Ui*U,G,V,Ord,Mod);
                    905:                }
                    906:                print(J,2);
                    907:                Int0 = Int1;
                    908:                S = cons(Ui,S);
                    909:        }
                    910:        print("");
                    911:        dp_gr_print(Dp);
                    912:        return reverse(S);
                    913: }
                    914:
                    915: def find_ssi2(C,G,Q,Rad,V,Ord) {
                    916:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                    917:        if ( type(Reduce=getopt(red)) == -1 ) Reduce = 0;
                    918:        for ( T = C, S = []; T != []; T = cdr(T) )
                    919:                if ( gen_nf(car(T),Q,V,Ord,Mod) ) S = cons(car(T),S);
                    920:        if ( S == [] )
                    921:                error("find_ssi2 : cannot happen");
                    922:        G1 = append(S,G);
                    923:        Int = ideal_intersection(G1,Q,V,Ord|mod=Mod);
                    924:        /* check whether (Q cap (G+S)) = G */
                    925:        if ( gen_gb_comp(Int,G,Mod) ) { print([0]); return reverse(S); }
                    926:
                    927: #if 0
                    928:        dp_ord(Ord); DC = map(dp_ptod,C,V);
1.9       ohara     929:        DC = qsort(DC,noro_pd.comp_tord); C = map(dp_dtop,DC,V);
1.1       noro      930: #else
1.9       ohara     931:        C = qsort(C,noro_pd.comp_tdeg);
1.1       noro      932: #endif
                    933:        if ( Reduce ) {
                    934:                for ( T = C, C1 = [], R1 = Rad; T != []; T = cdr(T) ) {
                    935:                        if ( !gen_nf(car(T),Rad,V,Ord,Mod) ) continue;
1.4       noro      936:                        if ( radical_membership(car(T),R1,V|mod=Mod)  ) {
1.1       noro      937:                                C1 = cons(car(T),C1);
                    938:                                R1 = cons(sq(car(T),Mod),R1);
                    939:                        }
                    940:                }
                    941:                print(["C",length(C),"->",length(C1)]);
                    942:                C = reverse(C1);
                    943:        }
                    944:        for ( T = C, S = []; T != []; T = cdr(T) ) {
                    945:                if ( !gen_nf(car(T),Rad,V,Ord,Mod) ) continue;
                    946:                Ui = U = car(T);
                    947:                S = cons([Ui,U],S);
                    948:        }
1.9       ohara     949:        S = qsort(S,noro_pd.comp_tdeg_first);
1.1       noro      950:        print("");
                    951:
                    952:        Dp = dp_gr_print(); dp_gr_print(0);
                    953:        Tmp = ttttt; TV = cons(Tmp,V); Ord1 = [[0,1],[Ord,length(V)]];
                    954:        Int0 = incremental_gb(append(vtol(ltov(G)*Tmp),vtol(ltov(Q)*(1-Tmp))),
                    955:                TV,Ord1|gbblock=[[0,length(G)]],mod=Mod);
                    956:        OK = [];
                    957:        while ( S != [] ) {
                    958:                Len = length(S); print("remaining gens : ",0); print(Len);
                    959:                S1 = [];
                    960:                for ( Start = Prev = 0; Start < Len; Start = Prev ) {
                    961:                        Cur = Start+1;
                    962:                        print(Start,2);
                    963:                        while ( Prev < Len ) {
                    964:                                for ( St = [], I = Prev; I < Cur; I++ ) St = cons(Tmp*S[I][0],St);
                    965:                                Int1 = nd_gr(append(Int0,St),TV,Mod,Ord1|gbblock=[[0,length(Int0)]],newelim=1);
                    966:                                if ( !Int1 ) {
                    967:                                        print("x",0); break;
                    968:                                }
                    969:                                Int = elimination(Int1,V);
                    970:                                if ( gen_gb_comp(Int,G,Mod) ) {
                    971:                                        print([Prev,Cur-1],2);
                    972:                                        Prev = Cur;
                    973:                                        Cur += (Prev-Start)+1;
                    974:                                        if ( Cur > Len ) Cur = Len;
                    975:                                        Int0 = Int1;
                    976:                                } else
                    977:                                        break;
                    978:                        }
                    979:                        for ( I = Start; I < Prev; I++ ) OK = cons(S[I][0],OK);
                    980:                        if ( Prev == Start ) {
                    981:                                Ui = S[I][0]; U = S[I][1];
                    982:                                Ui = gen_nf(Ui*U,G,V,Ord,Mod);
                    983:                                S1 = cons([Ui,U],S1);
                    984:                                Prev++;
                    985:                        }
                    986:                }
                    987:                S = reverse(S1);
                    988:                print("");
                    989:        }
                    990:        print("");
                    991:        OK = reverse(OK);
                    992:        dp_gr_print(Dp);
                    993:        return OK;
                    994: }
                    995:
                    996: /* SY primary decompsition */
                    997:
                    998: def sy_dec(B,V)
                    999: {
                   1000:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   1001:        if ( type(Lexdec=getopt(lexdec)) == -1 ) Lexdec = 0;
                   1002:        Ord = 0;
                   1003:        G = fast_gb(B,V,Mod,Ord);
                   1004:        Q = [];
                   1005:        IntQ = [1];
                   1006:        Gt = G;
                   1007:        First = 1;
                   1008:        while ( 1 ) {
                   1009:                if ( type(Gt[0]) == 1 ) break;
                   1010:                Pt = prime_dec(Gt,V|indep=1,lexdec=Lexdec,mod=Mod);
                   1011:                L = pseudo_dec(Gt,Pt,V,Ord|mod=Mod);
                   1012:                Qt = L[0]; Rt = L[1]; St = L[2];
                   1013:                IntQt = ideal_list_intersection(map(first,Qt),V,Ord|mod=Mod);
                   1014:                if ( First ) {
                   1015:                        IntQ = IntQt;
                   1016:                        Qi = Qt;
                   1017:                        First = 0;
                   1018:                } else {
                   1019:                        IntQ = ideal_intersection(IntQ,IntQt,V,Ord|mod=Mod);
                   1020:                        Q = append(Qt,Q);
                   1021:                }
                   1022:                if ( gen_gb_comp(IntQ,G,Mod) ) break;
                   1023:                for ( T = Rt; T != []; T = cdr(T) ) {
                   1024:                        if ( type(car(T)[0]) == 1 ) continue;
                   1025:                        U = sy_dec(car(T),V|lexdec=Lexdec,mod=Mod);
                   1026:                        IntQ = ideal_list_intersection(cons(IntQ,map(first,U)),
                   1027:                                V,Ord|mod=Mod);
                   1028:                        Q = append(U,Q);
                   1029:                        if ( gen_gb_comp(IntQ,G,Mod) ) break;
                   1030:                }
                   1031:                Gt = fast_gb(append(Gt,St),V,Mod,Ord);
                   1032:        }
                   1033:        Q = qd_remove_redundant_comp(G,Qi,Q,V,Ord|mod=Mod);
                   1034:        return append(Qi,Q);
                   1035: }
                   1036:
                   1037: def pseudo_dec(G,L,V,Ord)
                   1038: {
                   1039:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   1040:        N = length(L);
                   1041:        S = vector(N);
                   1042:        Q = vector(N);
                   1043:        R = vector(N);
                   1044:        L0 = map(first,L);
                   1045:        for ( I = 0; I < N; I++ ) {
                   1046:                LI = setminus(L0,[L0[I]]);
                   1047:                PI = ideal_list_intersection(LI,V,Ord|mod=Mod);
1.9       ohara    1048:                PI = qsort(PI,noro_pd.comp_tdeg);
1.1       noro     1049:                for ( T = PI; T != []; T = cdr(T) )
                   1050:                        if ( gen_nf(car(T),L0[I],V,Ord,Mod) ) break;
                   1051:                if ( T == [] ) error("separator : cannot happen");
1.11    ! noro     1052:                SI = sat_ind(G,car(T),V|mod=Mod);
1.1       noro     1053:                QI = SI[0];
                   1054:                S[I] = car(T)^SI[1];
                   1055:                PV = L[I][1];
                   1056:                V0 = setminus(V,PV);
                   1057: #if 0
                   1058:                GI = fast_gb(QI,append(V0,PV),Mod,
                   1059:                        [[Ord,length(V0)],[Ord,length(PV)]]);
                   1060: #else
                   1061:                GI = fast_gb(QI,append(V0,PV),Mod,
                   1062:                        [[2,length(V0)],[Ord,length(PV)]]);
                   1063: #endif
                   1064:                LCFI = lcfactor(GI,V0,Ord,Mod);
                   1065:                for ( F = 1, T = LCFI, Gt = QI; T != []; T = cdr(T) ) {
1.11    ! noro     1066:                        St = sat_ind(Gt,T[0],V|mod=Mod);
1.1       noro     1067:                        Gt = St[0]; F *= T[0]^St[1];
                   1068:                }
                   1069:                Q[I] = [Gt,L0[I]];
                   1070:                R[I] = fast_gb(cons(F,QI),V,Mod,Ord);
                   1071:        }
                   1072:        return [vtol(Q),vtol(R),vtol(S)];
                   1073: }
                   1074:
                   1075: def iso_comp(G,L,V,Ord)
                   1076: {
                   1077:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   1078:        if ( type(IsGB=getopt(isgb)) == -1 ) IsGB = 0;
                   1079:        if ( type(Iso=getopt(iso)) == -1 ) Iso = 0;
                   1080:        if ( type(Para=getopt(para)) == -1 ) Para = 0;
                   1081:        if ( type(Q=getopt(intq)) == -1 ) Q = 0;
1.4       noro     1082:        if ( type(S=getopt(sep)) == -1 ) S = 0;
1.1       noro     1083:
1.4       noro     1084:        if ( !S ) S = separator(L,V|mod=Mod);
1.1       noro     1085:        N = length(L);
                   1086:        print("comps : ",2); print(N); print("",2);
                   1087:        if ( Para ) {
                   1088:                Task = [];
                   1089:                for ( I = 0; I < N; I++ ) {
                   1090:                        T = ["noro_pd.locsat",G,V,L[I],S[I],Mod,IsGB,Iso,Q];
                   1091:                        Task = cons(T,Task);
                   1092:                }
                   1093:                Task = reverse(Task);
                   1094:                R = para_exec(Para,Task);
                   1095:                return R;
                   1096:        } else {
                   1097:                for ( I = 0, R = []; I < N; I++ ) {
                   1098:                        QI = locsat(G,V,L[I],S[I],Mod,IsGB,Iso,Q);
                   1099:                        if ( type(QI[0][0])==1 )
                   1100:                                error("iso_comp : cannot happen");
                   1101:                        print(".",2);
                   1102:                        R = cons(QI,R);
                   1103:                }
                   1104:                print("");
                   1105:                return reverse(R);
                   1106:        }
                   1107: }
                   1108:
                   1109: def locsat(G,V,L,S,Mod,IsGB,Iso,Q)
                   1110: {
                   1111:        P = L[0]; PV = L[1]; V0 = setminus(V,PV);
                   1112:        if ( Iso==1 ) {
                   1113:                QI = sat(G,S,V|isgb=IsGB,mod=Mod);
                   1114:                GI = elim_gb(QI,V0,PV,Mod,[[0,length(V0)],[0,length(PV)]]);
1.8       noro     1115:                GI = nd_gr(contraction(GI,V0|mod=Mod,allv=V),V,Mod,0);
1.1       noro     1116:        } else if ( Iso==0 ) {
                   1117:                HI = elim_gb(G,V0,PV,Mod,[[0,length(V0)],[0,length(PV)]]);
1.8       noro     1118:                GI = nd_gr(contraction(HI,V0|mod=Mod,allv=V),V,Mod,0);
1.1       noro     1119:                GI = sat(GI,S,V|isgb=IsGB,mod=Mod);
                   1120:        } else if ( Iso==2 ) {
                   1121:                HI = elim_gb(G,V0,PV,Mod,[[0,length(V0)],[0,length(PV)]]);
                   1122:                TV = ttttt;
                   1123:                if ( Mod )
                   1124:                        GI = nd_gr(cons(TV*S-1,HI),cons(TV,V0),Mod,[[0,1],[0,length(V0)]]);
                   1125:                else
                   1126:                        GI = nd_gr_trace(append(HI,[TV*S-1]),cons(TV,V0),
                   1127:                                1,1,[[0,1],[0,length(V0)]]|gbblock=[[0,length(HI)]]);
                   1128:                GI = elimination(GI,V);
1.8       noro     1129:                GI = nd_gr(contraction(GI,V0|mod=Mod,allv=V),V,Mod,0);
1.4       noro     1130:        } else if ( Iso==3 ) {
                   1131:                GI = sat(G,S,V|isgb=IsGB,mod=Mod);
1.1       noro     1132:        }
                   1133:        if ( Q )
                   1134:                GI = ideal_intersection(Q,GI,V,0|mod=Mod);
                   1135:        return [GI,P,PV];
                   1136: }
                   1137:
                   1138: /* GTZ primary decompsition */
                   1139:
                   1140: def prima_dec(B,V)
                   1141: {
                   1142:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   1143:        if ( type(Ord=getopt(ord)) == -1 ) Ord = 0;
                   1144:        G0 = fast_gb(B,V,Mod,0);
                   1145:        G = fast_gb(G0,V,Mod,Ord);
                   1146:        IntP = [1];
                   1147:        QD = [];
                   1148:        while ( 1 ) {
                   1149:                if ( type(G[0])==1 || ideal_inclusion(IntP,G0,V,0|mod=Mod) )
                   1150:                        break;
                   1151:                W = maxindep(G,V,Ord); NP = length(W);
                   1152:                V0 = setminus(V,W); N = length(V0);
                   1153:                V1 = append(V0,W);
                   1154:                G1 = fast_gb(G,V1,Mod,[[Ord,N],[Ord,NP]]);
                   1155:                LCF = lcfactor(G1,V0,Ord,Mod);
                   1156:                L = zprimacomp(G,V0|mod=Mod);
                   1157:                F = 1;
                   1158:                for ( T = LCF, G2 = G; T != []; T = cdr(T) ) {
1.11    ! noro     1159:                        S = sat_ind(G2,T[0],V1|mod=Mod);
1.1       noro     1160:                        G2 = S[0]; F *= T[0]^S[1];
                   1161:                }
                   1162:                for ( T = L, QL = []; T != []; T = cdr(T) )
                   1163:                        QL = cons(car(T)[0],QL);
                   1164:                Int = ideal_list_intersection(QL,V,0|mod=Mod);
                   1165:                IntP = ideal_intersection(IntP,Int,V,0|mod=Mod);
                   1166:                QD = append(QD,L);
                   1167:                F = gen_nf(F,G,V,0,Mod);
                   1168:                G = fast_gb(cons(F,G),V,Mod,Ord);
                   1169:        }
                   1170:        QD = qd_remove_redundant_comp(G0,[],QD,V,0);
                   1171:        return QD;
                   1172: }
                   1173:
                   1174: /* SL prime decomposition */
                   1175:
                   1176: def prime_dec(B,V)
                   1177: {
                   1178:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   1179:        if ( type(Indep=getopt(indep)) == -1 ) Indep = 0;
1.7       noro     1180:        if ( type(LexDec=getopt(lexdec)) == -1 ) LexDec = 0;
1.1       noro     1181:        if ( type(Rad=getopt(radical)) == -1 ) Rad = 0;
                   1182:        B = map(sq,B,Mod);
                   1183:        if ( LexDec )
                   1184:                PD = lex_predec1(B,V|mod=Mod);
                   1185:        else
                   1186:                PD = [B];
                   1187:        if ( length(PD) > 1 ) {
                   1188:                G = ideal_list_intersection(PD,V,0|mod=Mod);
                   1189:                PD = pd_remove_redundant_comp(G,PD,V,0|mod=Mod);
                   1190:        }
1.8       noro     1191:        R = []; RL = [];
                   1192:        for ( T = PD; T != []; T = cdr(T) ) {
                   1193:                PDT = prime_dec_main(car(T),V|indep=Indep,mod=Mod);
                   1194:                R = append(R,PDT[0]);
                   1195:                GT = nd_gr(PDT[1],V,Mod,0);
                   1196:                RL = append(RL,[GT]);
                   1197:        }
                   1198:        if ( LexDec ) R = pd_simp_comp(R,V|first=Indep,mod=Mod);
                   1199:        if ( Rad ) {
                   1200:                G = ideal_list_intersection(RL,V,0|mod=Mod);
                   1201:                return [R,G];
                   1202:        } else return R;
1.7       noro     1203: }
                   1204:
                   1205: def prime_dec2(B,V)
                   1206: {
                   1207:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   1208:        if ( type(Indep=getopt(indep)) == -1 ) Indep = 0;
                   1209:        if ( type(LexDec=getopt(lexdec)) == -1 ) LexDec = 0;
                   1210:        if ( type(Rad=getopt(radical)) == -1 ) Rad = 0;
                   1211:        if ( type(Para=getopt(para)) == -1 || type(Para) != 4 ) Para = [];
                   1212:        B = map(sq,B,Mod);
                   1213:        if ( LexDec )
                   1214:                PD = lex_predec1(B,V|mod=Mod);
                   1215:        else
                   1216:                PD = [B];
                   1217:        if ( length(PD) > 1 ) {
                   1218:                G = ideal_list_intersection(PD,V,0|mod=Mod);
                   1219:                PD = pd_remove_redundant_comp(G,PD,V,0|mod=Mod);
                   1220:        }
                   1221:        R = [];
                   1222:        for ( T = PD; T != []; T = cdr(T) )
                   1223:                R = append(prime_dec_main2(car(T),V|indep=Indep,mod=Mod,para=Para),R);
                   1224:        if ( Indep ) {
                   1225:                G = ideal_list_intersection(map(first,R),V,0|mod=Mod);
                   1226:                R = pd_simp_comp(R,V|first=1,mod=Mod);
                   1227:        } else {
                   1228:                G = ideal_list_intersection(R,V,0|mod=Mod);
                   1229:                R = pd_simp_comp(R,V|mod=Mod);
1.1       noro     1230:        }
                   1231:        return Rad ? [R,G] : R;
                   1232: }
                   1233:
1.8       noro     1234: /* returns [PD,rad(I)] */
                   1235:
1.1       noro     1236: def prime_dec_main(B,V)
                   1237: {
1.8       noro     1238:        Tpint = RTpint = 0;
1.1       noro     1239:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   1240:        if ( type(Indep=getopt(indep)) == -1 ) Indep = 0;
                   1241:        G = fast_gb(B,V,Mod,0);
                   1242:        IntP = [1];
                   1243:        PD = [];
1.4       noro     1244:        DG = ltov(map(dp_ptod,G,V));
                   1245:        for ( Ind = [], I = length(G)-1; I >= 0; I-- ) Ind = cons(I,Ind);
                   1246:        if ( Mod ) DG = map(dp_mod,DG,Mod,[]);
1.1       noro     1247:        while ( 1 ) {
1.8       noro     1248:                print([length(PD)],2);
1.1       noro     1249:                /* rad(G) subset IntP */
                   1250:                /* check if IntP subset rad(G) */
1.4       noro     1251:                /* print([length(PD),length(IntP)],2); */
1.7       noro     1252:                for ( T = IntP; T != []; T = cdr(T) )
1.8       noro     1253:                        if ( (G0 = radical_membership_sat(car(T),G,V|mod=Mod,isgb=1,dg=[DG,Ind])) ) {
1.1       noro     1254:                                F = car(T);
                   1255:                                break;
                   1256:                        }
1.8       noro     1257:                if ( T == [] ) {
                   1258:                        print(["pint",Tpint,"rpint",RTpint]);
                   1259:                        return [PD,IntP];
                   1260:                }
1.1       noro     1261:                PD0 = zprimecomp(G0,V,Indep|mod=Mod);
1.7       noro     1262:                Int = ideal_list_intersection(Indep?map(first,PD0):PD0,V,0|mod=Mod);
                   1263:                PD = append(PD,PD0);
1.8       noro     1264: #if 1
                   1265: T0=time();
1.4       noro     1266:                IntP = ideal_intersection_m(IntP,Int,V,0|mod=Mod);
1.8       noro     1267:                dp_ord(0); DC = map(dp_ptod,IntP,V);
1.9       ohara    1268:                DC = qsort(DC,noro_pd.comp_tord); IntP = map(dp_dtop,DC,V);
1.8       noro     1269: ACCUM_TIME(Tpint,RTpint)
1.4       noro     1270: #else
1.7       noro     1271:                IntP = ideal_intersection(IntP,Int,V,0|mod=Mod,gbblock=[[0,length(IntP)]]);
1.4       noro     1272: #endif
1.1       noro     1273:        }
                   1274: }
                   1275:
1.8       noro     1276: localf callsat,callzcomp;
                   1277:
                   1278: def callsat(F,G,V,Mod,DG)
                   1279: {
                   1280:        return radical_membership(F,G,V|mod=Mod,isgb=1,dg=DG,sat=1);
                   1281: }
                   1282:
                   1283: def callzcomp(F,V,Indep,Mod)
                   1284: {
                   1285:        PD0 = zprimecomp(F,V,Indep|mod=Mod);
                   1286:        Int = ideal_list_intersection(Indep?map(first,PD0):PD0,V,0|mod=Mod);
                   1287:        return [PD0,Int];
                   1288: }
                   1289:
1.7       noro     1290: def prime_dec_main2(B,V)
                   1291: {
                   1292:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   1293:        if ( type(Indep=getopt(indep)) == -1 ) Indep = 0;
                   1294:        if ( type(Para=getopt(para)) == -1 || type(Para) != 4 ) Para = [];
                   1295:        NPara = length(Para);
                   1296:
                   1297:        G = fast_gb(B,V,Mod,0);
                   1298:        IntP = [1];
                   1299:        PD = [];
                   1300:        DG = ltov(map(dp_ptod,G,V));
                   1301:        for ( Ind = [], I = length(G)-1; I >= 0; I-- ) Ind = cons(I,Ind);
                   1302:        if ( Mod ) DG = map(dp_mod,DG,Mod,[]);
1.8       noro     1303:        if ( NPara )
                   1304:                while ( 1 ) {
                   1305:                        IntPM = mingen(IntP,V);
                   1306:                        for ( T = IntPM, CallSat = []; T != []; T = cdr(T) )
                   1307:                                CallSat = cons(["noro_pd.callsat",car(T),G,V,Mod,[DG,Ind]],CallSat);
                   1308:                        CallSat = reverse(CallSat);
                   1309:                        /* SatL = [[..],0,[...],...] */
                   1310:                        SatL = para_exec(Para,CallSat);
                   1311:                        for ( T = SatL, Sat = []; T != []; T = cdr(T) ) if ( car(T) ) Sat = cons(car(T),Sat);
                   1312:                        if ( Sat == [] ) return PD;
                   1313:                        print(length(Sat),2); print("->",2);
                   1314:                        Sat = remove_identical_comp(Sat|mod=Mod);
                   1315:                        print(length(Sat));
                   1316:                        for ( T = Sat, CallComp = []; T != []; T = cdr(T) )
                   1317:                                CallComp = cons(["noro_pd.callzcomp",car(T),V,Indep,Mod],CallComp);
                   1318:                        CallComp = reverse(CallComp);
                   1319:                        /* PDL = [[PD0,Int],...] */
                   1320:                        PDL = para_exec(Para,CallComp);
                   1321:                        for ( T = PDL; T != []; T = cdr(T) ) PD = append(PD,car(T)[0]);
                   1322:                        Int = ideal_list_intersection(map(second,PDL),V,0|mod=Mod);
1.7       noro     1323:                        IntP = ideal_intersection(IntP,Int,V,0|mod=Mod,gbblock=[[0,length(IntP)]]);
                   1324:                }
1.8       noro     1325:        else
                   1326:                while ( 1 ) {
                   1327:                        /* rad(G) subset IntP */
                   1328:                        /* check if IntP subset rad(G) */
                   1329:                        /* print([length(PD),length(IntP)],2); */
                   1330:                        Sat = [];
                   1331:                        IntPM = mingen(IntP,V);
                   1332:                        for ( T = IntPM; T != [] && length(Sat) < 16; T = cdr(T) )
                   1333:                                if ( G0 = radical_membership(car(T),G,V|mod=Mod,isgb=1,dg=[DG,Ind],sat=1) )
                   1334:                                        Sat = cons(G0,Sat);
                   1335:                        if ( Sat == [] ) return PD;
                   1336:                        print(length(Sat),2); print("->",2);
                   1337:                        Sat = remove_identical_comp(Sat|mod=Mod);
                   1338:                        print(length(Sat));
                   1339:                        for ( T = Sat; T != []; T = cdr(T) ) {
                   1340:                                PD0 = zprimecomp(car(T),V,Indep|mod=Mod); PD = append(PD,PD0);
                   1341:                                Int = ideal_list_intersection(Indep?map(first,PD0):PD0,V,0|mod=Mod);
                   1342:                                IntP = ideal_intersection(IntP,Int,V,0|mod=Mod,gbblock=[[0,length(IntP)]]);
                   1343:                        }
                   1344:                }
1.7       noro     1345: }
                   1346:
1.1       noro     1347: /* pre-decomposition */
                   1348:
                   1349: def lex_predec1(B,V)
                   1350: {
                   1351:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   1352:        G = fast_gb(B,V,Mod,2);
                   1353:        for ( T = G; T != []; T = cdr(T) ) {
                   1354:                F = gen_fctr(car(T),Mod);
                   1355:                if ( length(F) > 2 || length(F) == 2 && F[1][1] > 1 ) {
                   1356:                        for ( R = [], S = cdr(F); S != []; S = cdr(S) ) {
                   1357:                                Ft = car(S)[0];
                   1358:                                Gt = map(ptozp,map(gen_nf,G,[Ft],V,0,Mod));
                   1359:                                R1 = fast_gb(cons(Ft,Gt),V,Mod,0);
                   1360:                                R = cons(R1,R);
                   1361:                        }
                   1362:                        return R;
                   1363:                }
                   1364:        }
                   1365:        return [G];
                   1366: }
                   1367:
                   1368: /* zero-dimensional prime/primary decomosition */
                   1369:
                   1370: def zprimedec(B,V,Mod)
                   1371: {
                   1372:        L = partial_decomp(B,V,Mod);
                   1373:        P = L[0]; NP = L[1];
                   1374:        R = [];
                   1375:        for ( ; P != []; P = cdr(P) ) R = cons(car(car(P)),R);
                   1376:        for ( T = NP; T != []; T = cdr(T) ) {
                   1377:                R1 = complete_decomp(car(T),V,Mod);
                   1378:                R = append(R1,R);
                   1379:        }
                   1380:        return R;
                   1381: }
                   1382:
                   1383: def zprimadec(B,V,Mod)
                   1384: {
                   1385:        L = partial_qdecomp(B,V,Mod);
                   1386:        Q = L[0]; NQ = L[1];
                   1387:        R = [];
                   1388:        for ( ; Q != []; Q = cdr(Q) ) {
                   1389:                T = car(Q); R = cons([T[0],T[1]],R);
                   1390:        }
                   1391:        for ( T = NQ; T != []; T = cdr(T) ) {
                   1392:                R1 = complete_qdecomp(car(T),V,Mod);
                   1393:                R = append(R1,R);
                   1394:        }
                   1395:        return R;
                   1396: }
                   1397:
                   1398: def complete_qdecomp(GD,V,Mod)
                   1399: {
                   1400:        GQ = GD[0]; GP = GD[1]; D = GD[2];
                   1401:        W = vars(GP);
                   1402:        PV = setminus(W,V);
                   1403:        N = length(V); PN = length(PV);
                   1404:        U = find_npos([GP,D],V,PV,Mod);
                   1405:        NV = ttttt;
                   1406:        M = gen_minipoly(cons(NV-U,GQ),cons(NV,V),PV,0,NV,Mod);
                   1407:        M = ppart(M,NV,Mod);
1.6       noro     1408:        MF = Mod ? modfctr(M,Mod) : fctr(M);
1.1       noro     1409:        R = [];
                   1410:        for ( T = cdr(MF); T != []; T = cdr(T) ) {
                   1411:                S = car(T);
                   1412:                Mt = subst(S[0],NV,U);
                   1413:                GP1 = fast_gb(cons(Mt,GP),W,Mod,0);
                   1414:                GQ1 = fast_gb(cons(Mt^S[1],GQ),W,Mod,0);
                   1415:                if ( PV != [] ) {
                   1416:                        GP1 = elim_gb(GP1,V,PV,Mod,[[0,N],[0,PN]]);
                   1417:                        GQ1 = elim_gb(GQ1,V,PV,Mod,[[0,N],[0,PN]]);
                   1418:                }
                   1419:                R = cons([GQ1,GP1],R);
                   1420:        }
                   1421:        return R;
                   1422: }
                   1423:
                   1424: def partial_qdecomp(B,V,Mod)
                   1425: {
                   1426:        Elim = (Elim=getopt(elim))&&type(Elim)!=-1 ? 1 : 0;
                   1427:        N = length(V);
                   1428:        W = vars(B);
                   1429:        PV = setminus(W,V);
                   1430:        NP = length(PV);
                   1431:        W = append(V,PV);
                   1432:        if ( Elim && PV != [] ) Ord = [[0,N],[0,NP]];
                   1433:        else Ord = 0;
                   1434:        if ( Mod )
                   1435:                B = nd_f4(B,W,Mod,Ord);
                   1436:        else
                   1437:                B = nd_gr_trace(B,W,1,GBCheck,Ord);
                   1438:        Q = []; NQ = [[B,B,vector(N+1)]];
                   1439:        for ( I = length(V)-1; I >= 0; I-- ) {
                   1440:                NQ1 = [];
                   1441:                for ( T = NQ; T != []; T = cdr(T) ) {
                   1442:                        L = partial_qdecomp0(car(T),V,PV,Ord,I,Mod);
                   1443:                        Q = append(L[0],Q);
                   1444:                        NQ1 = append(L[1],NQ1);
                   1445:                }
                   1446:                NQ = NQ1;
                   1447:        }
                   1448:        return [Q,NQ];
                   1449: }
                   1450:
                   1451: def partial_qdecomp0(GD,V,PV,Ord,I,Mod)
                   1452: {
                   1453:        GQ = GD[0]; GP = GD[1]; D = GD[2];
                   1454:        N = length(V); PN = length(PV);
                   1455:        W = append(V,PV);
                   1456:        VI = V[I];
                   1457:        M = gen_minipoly(GQ,V,PV,Ord,VI,Mod);
                   1458:        M = ppart(M,VI,Mod);
                   1459:        if ( Mod )
                   1460:                MF = modfctr(M,Mod);
                   1461:        else
                   1462:                MF = fctr(M);
                   1463:        Q = []; NQ = [];
                   1464:        if ( length(MF) == 2 && MF[1][1] == 1 ) {
                   1465:                D1 = D*1; D1[I] = M;
                   1466:                GQelim = elim_gb(GQ,V,PV,Mod,Ord);
                   1467:                GPelim = elim_gb(GP,V,PV,Mod,Ord);
                   1468:                LD = ldim(GQelim,V);
                   1469:                if ( deg(M,VI) == LD )
                   1470:                        Q = cons([GQelim,GPelim,D1],Q);
                   1471:                else
                   1472:                        NQ = cons([GQelim,GPelim,D1],NQ);
                   1473:                return [Q,NQ];
                   1474:        }
                   1475:        for ( T = cdr(MF); T != []; T = cdr(T) ) {
                   1476:                S = car(T); Mt = S[0]; D1 = D*1; D1[I] = Mt;
                   1477:
                   1478:                GQ1 = fast_gb(cons(Mt^S[1],GQ),W,Mod,Ord);
                   1479:                GQelim = elim_gb(GQ1,V,PV,Mod,Ord);
                   1480:                GP1 = fast_gb(cons(Mt,GP),W,Mod,Ord);
                   1481:                GPelim = elim_gb(GP1,V,PV,Mod,Ord);
                   1482:
                   1483:                D1[N] = LD = ldim(GPelim,V);
                   1484:
                   1485:                for ( J = 0; J < N; J++ )
                   1486:                        if ( D1[J] && deg(D1[J],V[J]) == LD ) break;
                   1487:                if ( J < N )
                   1488:                        Q = cons([GQelim,GPelim,D1],Q);
                   1489:                else
                   1490:                        NQ = cons([GQelim,GPelim,D1],NQ);
                   1491:        }
                   1492:        return [Q,NQ];
                   1493: }
                   1494:
                   1495: def complete_decomp(GD,V,Mod)
                   1496: {
                   1497:        G = GD[0]; D = GD[1];
                   1498:        W = vars(G);
                   1499:        PV = setminus(W,V);
                   1500:        N = length(V); PN = length(PV);
                   1501:        U = find_npos(GD,V,PV,Mod);
                   1502:        NV = ttttt;
                   1503:        M = gen_minipoly(cons(NV-U,G),cons(NV,V),PV,0,NV,Mod);
                   1504:        M = ppart(M,NV,Mod);
1.6       noro     1505:        MF = Mod ? modfctr(M,Mod) : fctr(M);
1.1       noro     1506:        if ( length(MF) == 2 ) return [G];
                   1507:        R = [];
                   1508:        for ( T = cdr(MF); T != []; T = cdr(T) ) {
                   1509:                Mt = subst(car(car(T)),NV,U);
                   1510:                G1 = fast_gb(cons(Mt,G),W,Mod,0);
                   1511:                if ( PV != [] ) G1 = elim_gb(G1,V,PV,Mod,[[0,N],[0,PN]]);
                   1512:                R = cons(G1,R);
                   1513:        }
                   1514:        return R;
                   1515: }
                   1516:
                   1517: def partial_decomp(B,V,Mod)
                   1518: {
                   1519:        Elim = (Elim=getopt(elim))&&type(Elim)!=-1 ? 1 : 0;
                   1520:        N = length(V);
                   1521:        W = vars(B);
                   1522:        PV = setminus(W,V);
                   1523:        NP = length(PV);
                   1524:        W = append(V,PV);
                   1525:        if ( Elim && PV != [] ) Ord = [[0,N],[0,NP]];
                   1526:        else Ord = 0;
                   1527:        if ( Mod )
                   1528:                B = nd_f4(B,W,Mod,Ord);
                   1529:        else
                   1530:                B = nd_gr_trace(B,W,1,GBCheck,Ord);
                   1531:        P = []; NP = [[B,vector(N+1)]];
                   1532:        for ( I = length(V)-1; I >= 0; I-- ) {
                   1533:                NP1 = [];
                   1534:                for ( T = NP; T != []; T = cdr(T) ) {
                   1535:                        L = partial_decomp0(car(T),V,PV,Ord,I,Mod);
                   1536:                        P = append(L[0],P);
                   1537:                        NP1 = append(L[1],NP1);
                   1538:                }
                   1539:                NP = NP1;
                   1540:        }
                   1541:        return [P,NP];
                   1542: }
                   1543:
                   1544: def partial_decomp0(GD,V,PV,Ord,I,Mod)
                   1545: {
                   1546:        G = GD[0]; D = GD[1];
                   1547:        N = length(V); PN = length(PV);
                   1548:        W = append(V,PV);
                   1549:        VI = V[I];
                   1550:        M = gen_minipoly(G,V,PV,Ord,VI,Mod);
                   1551:        M = ppart(M,VI,Mod);
                   1552:        if ( Mod )
                   1553:                MF = modfctr(M,Mod);
                   1554:        else
                   1555:                MF = fctr(M);
                   1556:        if ( length(MF) == 2 && MF[1][1] == 1 ) {
                   1557:                D1 = D*1;
                   1558:                D1[I] = M;
                   1559:                Gelim = elim_gb(G,V,PV,Mod,Ord);
                   1560:                D1[N] = LD = ldim(Gelim,V);
                   1561:                GD1 = [Gelim,D1];
                   1562:                for ( J = 0; J < N; J++ )
                   1563:                        if ( D1[J] && deg(D1[J],V[J]) == LD )
                   1564:                                return [[GD1],[]];
                   1565:                return [[],[GD1]];
                   1566:        }
                   1567:        P = []; NP = [];
                   1568:        GI = elim_gb(G,V,PV,Mod,Ord);
                   1569:        for ( T = cdr(MF); T != []; T = cdr(T) ) {
                   1570:                Mt = car(car(T));
                   1571:                D1 = D*1;
                   1572:                D1[I] = Mt;
                   1573:                GIt = map(gen_nf,GI,[Mt],V,Ord,Mod);
                   1574:                G1 = cons(Mt,GIt);
                   1575:                Gelim = elim_gb(G1,V,PV,Mod,Ord);
                   1576:                D1[N] = LD = ldim(Gelim,V);
                   1577:                for ( J = 0; J < N; J++ )
                   1578:                        if ( D1[J] && deg(D1[J],V[J]) == LD ) break;
                   1579:                if ( J < N )
                   1580:                        P = cons([Gelim,D1],P);
                   1581:                else
                   1582:                        NP = cons([Gelim,D1],NP);
                   1583:        }
                   1584:        return [P,NP];
                   1585: }
                   1586:
                   1587: /* prime/primary components over rational function field */
                   1588:
                   1589: def zprimacomp(G,V) {
                   1590:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   1591:        L = zprimadec(G,V,0|mod=Mod);
                   1592:        R = [];
                   1593:        dp_ord(0);
                   1594:        for ( T = L; T != []; T = cdr(T) ) {
                   1595:                S = car(T);
                   1596:                UQ = contraction(S[0],V|mod=Mod);
                   1597:                UP = contraction(S[1],V|mod=Mod);
                   1598:                R = cons([UQ,UP],R);
                   1599:        }
                   1600:        return R;
                   1601: }
                   1602:
                   1603: def zprimecomp(G,V,Indep) {
                   1604:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   1605:        W = maxindep(G,V,0|mod=Mod);
                   1606:        V0 = setminus(V,W);
                   1607:        V1 = append(V0,W);
                   1608: #if 0
                   1609:        O1 = [[0,length(V0)],[0,length(W)]];
                   1610:        G1 = fast_gb(G,V1,Mod,O1);
                   1611:        dp_ord(0);
                   1612: #else
                   1613:        G1 = G;
                   1614: #endif
                   1615:        PD = zprimedec(G1,V0,Mod);
                   1616:        dp_ord(0);
                   1617:        R = [];
                   1618:        for ( T = PD; T != []; T = cdr(T) ) {
                   1619:                U = contraction(car(T),V0|mod=Mod);
                   1620:                U = nd_gr(U,V,Mod,0);
1.7       noro     1621:                R = cons(Indep?[U,W]:U,R);
1.1       noro     1622:        }
1.7       noro     1623:        return R;
1.1       noro     1624: }
                   1625:
                   1626: def fast_gb(B,V,Mod,Ord)
                   1627: {
                   1628:        if ( type(Block=getopt(gbblock)) == -1 ) Block = 0;
                   1629:        if ( type(NoRA=getopt(nora)) == -1 ) NoRA = 0;
                   1630:        if ( type(Trace=getopt(trace)) == -1 ) Trace = 0;
                   1631:        if ( Mod )
                   1632:                G = nd_f4(B,V,Mod,Ord|nora=NoRA);
                   1633:        else if ( F4 )
                   1634:                G = map(ptozp,f4_chrem(B,V,Ord));
                   1635:        else if ( Trace ) {
                   1636:                if ( Block )
                   1637:                        G = nd_gr_trace(B,V,1,1,Ord|nora=NoRA,gbblock=Block);
                   1638:                else
                   1639:                        G = nd_gr_trace(B,V,1,1,Ord|nora=NoRA);
                   1640:        } else {
                   1641:                if ( Block )
                   1642:                        G = nd_gr(B,V,0,Ord|nora=NoRA,gbblock=Block);
                   1643:                else
                   1644:                        G = nd_gr(B,V,0,Ord|nora=NoRA);
                   1645:        }
                   1646:        return G;
                   1647: }
                   1648:
                   1649: def incremental_gb(A,V,Ord)
                   1650: {
                   1651:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   1652:        if ( type(Block=getopt(gbblock)) == -1 ) Block = 0;
                   1653:        if ( Mod ) {
                   1654:                if ( Block )
                   1655:                        G = nd_gr(A,V,Mod,Ord|gbblock=Block);
                   1656:                else
                   1657:                        G = nd_gr(A,V,Mod,Ord);
                   1658:        } else if ( Procs ) {
                   1659:                Arg0 = ["nd_gr",A,V,0,Ord];
                   1660:                Arg1 = ["nd_gr_trace",A,V,1,GBCheck,Ord];
                   1661:                G = competitive_exec(Procs,Arg0,Arg1);
                   1662:        } else if ( Block )
                   1663:                G = nd_gr(A,V,0,Ord|gbblock=Block);
                   1664:        else
                   1665:                G = nd_gr(A,V,0,Ord);
                   1666:        return G;
                   1667: }
                   1668:
                   1669: def elim_gb(G,V,PV,Mod,Ord)
                   1670: {
                   1671:        N = length(V); PN = length(PV);
                   1672:        O1 = [[0,N],[0,PN]];
                   1673:        if ( Ord == O1 )
                   1674:                Ord = Ord[0][0];
                   1675:        if ( Mod ) /* XXX */ {
                   1676:                for ( T = G, H = []; T != []; T = cdr(T) )
                   1677:                        if ( car(T) ) H = cons(car(T),H);
                   1678:                G = reverse(H);
                   1679:                G = dp_gr_mod_main(G,V,0,Mod,Ord);
                   1680:        } else if ( EProcs ) {
                   1681: #if 1
                   1682:                Arg0 = ["dp_gr_main",G,V,0,0,Ord];
                   1683: #else
                   1684:                Arg0 = ["nd_gr",G,V,0,Ord];
                   1685: #endif
                   1686:                Arg1 = ["noro_pd.nd_gr_rat",G,V,PV,O1,Ord];
                   1687:                G = competitive_exec(EProcs,Arg0,Arg1);
                   1688:        } else if ( GBRat ) {
                   1689:                G1 = nd_gr(G,append(V,PV),0,O1);
1.8       noro     1690:                if ( GBRat == 1 )
                   1691:                G1 = nd_gr_postproc(G1,V,0,Ord,0|nora=1);
1.1       noro     1692:                return G1;
                   1693:        } else
                   1694: #if 1
1.2       noro     1695: #if 0
1.1       noro     1696:                G = dp_gr_main(G,V,0,0,Ord);
                   1697: #else
                   1698:                G = nd_gr_trace(G,V,1,1,Ord);
                   1699: #endif
                   1700: #else
                   1701:                G = nd_gr(G,V,0,Ord);
                   1702: #endif
                   1703:        return G;
                   1704: }
                   1705:
                   1706: def ldim(G,V)
                   1707: {
                   1708:        O0 = dp_ord(); dp_ord(0);
                   1709:        D = length(dp_mbase(map(dp_ptod,G,V)));
                   1710:        dp_ord(O0);
                   1711:        return D;
                   1712: }
                   1713:
                   1714: /* over Q only */
                   1715:
                   1716: def make_mod_subst(GD,V,PV,HC)
                   1717: {
                   1718:        N = length(V);
                   1719:        PN = length(PV);
                   1720:        G = GD[0]; D = GD[1];
                   1721:        for ( I = 0; ; I = (I+1)%100 ) {
                   1722:                Mod = lprime(I);
                   1723:                S = [];
                   1724:                for ( J = PN-1; J >= 0; J-- )
                   1725:                        S = append([PV[J],random()%Mod],S);
                   1726:                for ( T = HC; T != []; T = cdr(T) )
                   1727:                        if ( !(subst(car(T),S)%Mod) ) break;
                   1728:                if ( T != [] ) continue;
                   1729:                for ( J = 0; J < N; J++ ) {
                   1730:                        M = subst(D[J],S);
                   1731:                        F = modsqfr(M,Mod);
                   1732:                        if ( length(F) != 2 || F[1][1] != 1 ) break;
                   1733:                }
                   1734:                if ( J < N ) continue;
                   1735:                G0 = map(subst,G,S);
                   1736:                return [G0,Mod];
                   1737:        }
                   1738: }
                   1739:
                   1740: def rsgn()
                   1741: {
                   1742:        return random()%2 ? 1 : -1;
                   1743: }
                   1744:
                   1745: def find_npos(GD,V,PV,Mod)
                   1746: {
                   1747:        N = length(V); PN = length(PV);
                   1748:        G = GD[0]; D = GD[1]; LD = D[N];
1.5       noro     1749:        DH = map(dp_dtop,map(dp_ht,map(dp_ptod,D,V)),V);
1.1       noro     1750:        Ord0 = dp_ord(); dp_ord(0);
                   1751:        HC = map(dp_hc,map(dp_ptod,G,V));
                   1752:        dp_ord(Ord0);
                   1753:        if ( !Mod ) {
                   1754:                W = append(V,PV);
                   1755:                G1 = nd_gr_trace(G,W,1,GBCheck,[[0,N],[0,PN]]);
                   1756:                L = make_mod_subst([G1,D],V,PV,HC);
                   1757:                return find_npos([L[0],D],V,[],L[1]);
                   1758:        }
                   1759:        N = length(V);
                   1760:        NV = ttttt;
                   1761:        for ( B = 2; ; B++ ) {
                   1762:                for ( J = N-2; J >= 0; J-- ) {
1.5       noro     1763:                        for ( U = 0, K = J; K < N; K++ ) {
                   1764:                                if ( DH[K] == V[K] ) continue;
1.1       noro     1765:                                U += rsgn()*((random()%B+1))*V[K];
1.5       noro     1766:                        }
1.6       noro     1767: #if 0
1.1       noro     1768:                        M = minipolym(G,V,0,U,NV,Mod);
1.6       noro     1769: #else
                   1770:                        M = gen_minipoly(cons(NV-U,G),cons(NV,V),PV,0,NV,Mod);
                   1771: #endif
1.1       noro     1772:                        if ( deg(M,NV) == LD ) return U;
                   1773:                }
                   1774:        }
                   1775: }
                   1776:
                   1777: def gen_minipoly(G,V,PV,Ord,VI,Mod)
                   1778: {
1.6       noro     1779:        O0 = dp_ord();
1.1       noro     1780:        if ( PV == [] ) {
                   1781:                NV = sssss;
                   1782:                if ( Mod )
                   1783:                        M = minipolym(G,V,Ord,VI,NV,Mod);
                   1784:                else
                   1785:                        M = minipoly(G,V,Ord,VI,NV);
1.6       noro     1786:                dp_ord(O0);
1.1       noro     1787:                return subst(M,NV,VI);
                   1788:        }
                   1789:        W = setminus(V,[VI]);
                   1790:        PV1 = cons(VI,PV);
                   1791: #if 0
                   1792:        while ( 1 ) {
                   1793:                V1 = append(W,PV1);
                   1794:                if ( Mod )
                   1795:                        G = nd_gr(G,V1,Mod,[[0,1],[0,length(V1)-1]]|nora=1);
                   1796:                else
                   1797:                        G = nd_gr_trace(G,V1,1,GBCheck,[[0,1],[0,length(V1)-1]]|nora=1);
                   1798:                if ( W == [] ) break;
                   1799:                else {
                   1800:                        W = cdr(W);
                   1801:                        G = elimination(G,cdr(V1));
                   1802:                }
                   1803:        }
1.8       noro     1804: #elif 1
1.1       noro     1805:        if ( Mod ) {
                   1806:                V1 = append(W,PV1);
                   1807:                G = nd_gr(G,V1,Mod,[[0,length(W)],[0,length(PV1)]]);
                   1808:                G = elimination(G,PV1);
                   1809:        } else {
                   1810:                PV2 = setminus(PV1,[PV1[length(PV1)-1]]);
                   1811:                V2 = append(W,PV2);
                   1812:                G = nd_gr_trace(G,V2,1,GBCheck,[[0,length(W)],[0,length(PV2)]]|nora=1);
                   1813:                G = elimination(G,PV1);
                   1814:        }
                   1815: #else
                   1816:        V1 = append(W,PV1);
                   1817:        if ( Mod )
                   1818:                G = nd_gr(G,V1,Mod,[[0,length(W)],[0,length(PV1)]]|nora=1);
                   1819:        else
                   1820:                G = nd_gr_trace(G,V1,1,GBCheck,[[0,length(W)],[0,length(PV1)]]|nora=1);
                   1821:        G = elimination(G,PV1);
                   1822: #endif
                   1823:        if ( Mod )
                   1824:                G = nd_gr(G,PV1,Mod,[[0,1],[0,length(PV)]]|nora=1);
                   1825:        else
                   1826:                G = nd_gr_trace(G,PV1,1,GBCheck,[[0,1],[0,length(PV)]]|nora=1);
                   1827:        for ( M = car(G), T = cdr(G); T != []; T = cdr(T) )
                   1828:                if ( deg(car(T),VI) < deg(M,VI) ) M = car(T);
1.6       noro     1829:        dp_ord(O0);
1.1       noro     1830:        return M;
                   1831: }
                   1832:
                   1833: def indepset(V,H)
                   1834: {
                   1835:        if ( H == [] ) return V;
                   1836:        N = -1;
                   1837:        for ( T = V; T != []; T = cdr(T) ) {
                   1838:                VI = car(T);
                   1839:                HI = [];
                   1840:                for ( S = H; S != []; S = cdr(S) )
                   1841:                        if ( !tdiv(car(S),VI) ) HI = cons(car(S),HI);
                   1842:                RI = indepset(setminus(V,[VI]),HI);
                   1843:                if ( length(RI) > N ) {
                   1844:                        R = RI; N = length(RI);
                   1845:                }
                   1846:        }
                   1847:        return R;
                   1848: }
                   1849:
                   1850: def maxindep(B,V,O)
                   1851: {
                   1852:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   1853:        G = fast_gb(B,V,Mod,O);
                   1854:        Old = dp_ord();
                   1855:        dp_ord(O);
                   1856:        H = map(dp_dtop,map(dp_ht,map(dp_ptod,G,V)),V);
                   1857:        H = map(sq,H,0);
                   1858:        H = nd_gr(H,V,0,0);
                   1859:        H = monodec0(H,V);
                   1860:        N = length(V);
                   1861:        Dep = [];
                   1862:        for ( T = H, Len = N+1; T != []; T = cdr(T) ) {
                   1863:                M = length(car(T));
                   1864:                if ( M < Len ) {
                   1865:                        Dep = [car(T)];
                   1866:                        Len = M;
                   1867:                } else if ( M == Len )
                   1868:                        Dep = cons(car(T),Dep);
                   1869:        }
                   1870:        R = setminus(V,Dep[0]);
                   1871:        dp_ord(Old);
                   1872:        return R;
                   1873: }
                   1874:
1.7       noro     1875: def maxindep2(B,V,O)
                   1876: {
                   1877:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   1878:        G = fast_gb(B,V,Mod,O);
                   1879:        Old = dp_ord();
                   1880:        dp_ord(O);
                   1881:        H = map(dp_dtop,map(dp_ht,map(dp_ptod,G,V)),V);
                   1882:        H = map(sq,H,0);
                   1883:        H = nd_gr(H,V,0,0);
                   1884:        H = monodec0(H,V);
                   1885:        N = length(V);
                   1886:        Dep = [];
                   1887:        for ( T = H, Len = N+1; T != []; T = cdr(T) ) {
                   1888:                M = length(car(T));
                   1889:                if ( M < Len ) {
                   1890:                        Dep = [car(T)];
                   1891:                        Len = M;
                   1892:                } else if ( M == Len )
                   1893:                        Dep = cons(car(T),Dep);
                   1894:        }
                   1895:        R = [];
                   1896:        for ( T = Dep; T != []; T = cdr(T) )
                   1897:                R = cons(setminus(V,car(T)),R);
                   1898:        dp_ord(Old);
                   1899:        return reverse(R);
                   1900: }
                   1901:
                   1902:
1.1       noro     1903: /* ideal operations */
                   1904: def contraction(G,V)
                   1905: {
1.8       noro     1906:        if ( type(AllV=getopt(allv)) == -1 ) AllV = 0;
1.1       noro     1907:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
1.8       noro     1908:
                   1909:        if ( RepColon ) return contraction_m(G,V|allv=AllV,mod=Mod);
                   1910:
1.1       noro     1911:        C = [];
                   1912:        for ( T = G; T != []; T = cdr(T) ) {
                   1913:                C1 = dp_hc(dp_ptod(car(T),V));
                   1914:                S = gen_fctr(C1,Mod);
                   1915:                for ( S = cdr(S); S != []; S = cdr(S) )
                   1916:                        if ( !member(S[0][0],C) ) C = cons(S[0][0],C);
                   1917:        }
                   1918:        W = vars(G);
                   1919:        PV = setminus(W,V);
1.8       noro     1920:        if ( AllV ) W = AllV;
                   1921:        else W = append(V,PV);
1.1       noro     1922:        NV = ttttt;
1.8       noro     1923:        if ( SuccSat ) {
                   1924:                W1 = cons(NV,W);
                   1925:                O1 = [[0,1],[0,length(W)]];
                   1926:                Block = [];
                   1927:                for ( T = C; T != []; T = cdr(T) ) {
                   1928:                        G1 = nd_gr(append(G,[NV*car(T)-1]),W1,Mod,O1|gbblock=Block);
                   1929:                        G = elimination(G1,W);
                   1930:                        Block = [[0,length(G)]];
                   1931:                }
                   1932:        } else {
                   1933:                for ( T = C, S = 1; T != []; T = cdr(T) )
                   1934:                        S *= car(T);
                   1935:                G = saturation([G,NV],S,W|mod=Mod);
                   1936:        }
1.1       noro     1937:        return G;
                   1938: }
                   1939:
1.8       noro     1940: def contraction_m(G,V)
                   1941: {
                   1942:        if ( type(AllV=getopt(allv)) == -1 ) AllV = 0;
                   1943:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   1944:        C = [];
                   1945:        for ( T = G; T != []; T = cdr(T) ) {
                   1946:                C1 = dp_hc(dp_ptod(car(T),V));
                   1947:                S = gen_fctr(C1,Mod);
                   1948:                for ( S = cdr(S); S != []; S = cdr(S) )
                   1949:                        if ( !member(S[0][0],C) ) C = cons(S[0][0],C);
                   1950:        }
                   1951:        W = vars(G);
                   1952:        PV = setminus(W,V);
                   1953:        if ( AllV ) W = AllV;
                   1954:        else W = append(V,PV);
                   1955:        H = H0 = G;
                   1956:        while ( 1 ) {
                   1957:                for ( T = C; T != []; T = cdr(T) )
                   1958:                        H = map(sdiv,ideal_intersection_m([car(T)],H,W,0),car(T));
                   1959:                H = nd_gr(H,W,0,0);
                   1960:                if ( gb_comp(H0,H) ) break;
                   1961:                else H0 = H;
                   1962:        }
                   1963:        return H;
                   1964: }
                   1965:
1.1       noro     1966: def ideal_list_intersection(L,V,Ord)
                   1967: {
                   1968:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
1.4       noro     1969:        if ( type(IsGB=getopt(isgb)) == -1 ) IsGB = 0;
1.1       noro     1970:        N = length(L);
                   1971:        if ( N == 0 ) return [1];
1.4       noro     1972:        if ( N == 1 )
                   1973:                return IsGB ? L[0] : fast_gb(L[0],V,Mod,Ord);
1.8       noro     1974:        else {
1.4       noro     1975:                for ( I = 0, T = [1]; I < N; I++ )
                   1976:                        T = ideal_intersection_m(T,L[I],V,Ord|mod=Mod);
                   1977:                T = nd_gr(T,V,Mod,Ord);
                   1978:                return T;
1.1       noro     1979:        }
                   1980: }
                   1981:
1.4       noro     1982: def call_ideal_list_intersection(L,V,Mod,Ord,IsGB)
1.1       noro     1983: {
1.4       noro     1984:        return ideal_list_intersection(L,V,Ord|mod=Mod,isgb=IsGB);
1.1       noro     1985: }
                   1986:
                   1987: def ideal_intersection(A,B,V,Ord)
                   1988: {
                   1989:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   1990:        if ( type(Block=getopt(gbblock)) == -1 ) Block = 0;
                   1991:        T = ttttt;
                   1992:        if ( Mod ) {
                   1993:                if ( Block )
                   1994:                        G = nd_gr(append(vtol(ltov(A)*T),vtol(ltov(B)*(1-T))),
                   1995:                                cons(T,V),Mod,[[0,1],[Ord,length(V)]]|gbblock=Block,nora=0);
                   1996:                else
                   1997:                        G = nd_gr(append(vtol(ltov(A)*T),vtol(ltov(B)*(1-T))),
                   1998:                                cons(T,V),Mod,[[0,1],[Ord,length(V)]]|nora=0);
                   1999:        } else
                   2000:        if ( Procs ) {
                   2001:                Arg0 = ["nd_gr",
                   2002:                        append(vtol(ltov(A)*T),vtol(ltov(B)*(1-T))),
                   2003:                        cons(T,V),0,[[0,1],[Ord,length(V)]]];
                   2004:                Arg1 = ["nd_gr_trace",
                   2005:                        append(vtol(ltov(A)*T),vtol(ltov(B)*(1-T))),
                   2006:                        cons(T,V),1,GBCheck,[[0,1],[Ord,length(V)]]];
                   2007:                G = competitive_exec(Procs,Arg0,Arg1);
                   2008:        } else {
                   2009:                if ( Block )
                   2010:                        G = nd_gr(append(vtol(ltov(A)*T),vtol(ltov(B)*(1-T))),
                   2011:                                cons(T,V),0,[[0,1],[Ord,length(V)]]|gbblock=Block,nora=0);
                   2012:                else
                   2013:                        G = nd_gr(append(vtol(ltov(A)*T),vtol(ltov(B)*(1-T))),
                   2014:                                cons(T,V),0,[[0,1],[Ord,length(V)]]|nora=0);
                   2015:        }
                   2016:        G0 = elimination(G,V);
                   2017:        if ( 0 && !Procs )
                   2018:                G0 = nd_gr_postproc(G0,V,Mod,Ord,0);
                   2019:        return G0;
                   2020: }
                   2021:
1.4       noro     2022:
                   2023: def aa(A) { return [A,A]; }
                   2024:
                   2025: def ideal_intersection_m(A,B,V,Ord)
                   2026: {
                   2027:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2028:
                   2029:        dp_ord(Ord);
                   2030:        DA = map(dp_ptod,A,V); DB = ltov(map(dp_ptod,B,V));
                   2031:        if ( Mod ) {
                   2032:                DA = map(dp_mod,DA,Mod,[]); DB = map(dp_mod,DB,Mod,[]);
                   2033:                setmod(Mod);
                   2034:        }
                   2035:        N = length(B);
                   2036:        for ( Ind = [], I = N-1; I >= 0; I-- ) Ind = cons(I,Ind);
                   2037:        for ( T = DA, C = []; T != []; T = cdr(T) ) {
                   2038:                L = Mod?dp_true_nf_mod(Ind,car(T),DB,1,Mod):dp_true_nf(Ind,car(T),DB,1);
                   2039:                R = dp_dtop(L[0],V); Q = dp_dtop(car(T)*L[1]-L[0],V);
                   2040:                C = cons([R,-Q],C);
                   2041:        }
                   2042:        G = nd_gr(append(C,map(aa,B)),V,Mod,[1,Ord]|intersect=1);
                   2043:        G = map(second,G);
                   2044:        return G;
                   2045: }
                   2046:
1.1       noro     2047: /* returns GB if F notin rad(G) */
                   2048:
                   2049: def radical_membership(F,G,V) {
                   2050:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2051:        if ( type(IsGB=getopt(isgb)) == -1 ) IsGB = 0;
1.4       noro     2052:        if ( type(L=getopt(dg)) == -1 ) L = 0;
1.7       noro     2053:        if ( type(Sat=getopt(sat)) == -1 ) Sat = 0;
1.4       noro     2054:        dp_ord(0);
                   2055:        if ( L ) { DG = L[0]; Ind = L[1]; }
                   2056:        else {
                   2057:                DG = ltov(map(dp_ptod,G,V));
                   2058:                if ( Mod ) DG = map(dp_mod,DG,Mod,[]);
                   2059:                for ( Ind = [], I = length(G)-1; I >= 0; I-- ) Ind = cons(I,Ind);
                   2060:        }
                   2061:        DF = dp_ptod(F,V); DFI = dp_ptod(1,V);
                   2062:        if ( Mod ) {
                   2063:                DF = dp_mod(DF,Mod,[]); DFI = dp_mod(DFI,Mod,[]);
                   2064:                setmod(Mod);
                   2065:        }
                   2066:        for ( I = 0; I < 3; I++ ) {
                   2067:                DFI = Mod?dp_nf_mod(Ind,DF*DFI,DG,0,Mod):dp_nf(Ind,DF*DFI,DG,0);
                   2068:                if ( !DFI ) return 0;
                   2069:        }
1.1       noro     2070:        NV = ttttt;
                   2071:        if ( IsGB )
                   2072:                T = nd_gr(append(G,[NV*F-1]),cons(NV,V),Mod,0
                   2073:                        |gbblock=[[0,length(G)]]);
                   2074:        else
                   2075:                T = nd_gr(append(G,[NV*F-1]),cons(NV,V),Mod,0);
1.7       noro     2076:        if ( type(car(T)) == 1 ) return 0;
                   2077:        else if ( Sat ) {
                   2078:                G1 = fast_gb(T,cons(NV,V),Mod,[[0,1],[0,length(V)]]);
                   2079:                G0 = elimination(G1,V);
                   2080:                return G0;
                   2081:        } else return [T,NV];
1.1       noro     2082: }
                   2083:
1.8       noro     2084: def radical_membership_sat(F,G,V) {
                   2085:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2086:        if ( type(IsGB=getopt(isgb)) == -1 ) IsGB = 0;
                   2087:        if ( type(L=getopt(dg)) == -1 ) L = 0;
                   2088:        dp_ord(0);
                   2089:        if ( L ) { DG = L[0]; Ind = L[1]; }
                   2090:        else {
                   2091:                DG = ltov(map(dp_ptod,G,V));
                   2092:                if ( Mod ) DG = map(dp_mod,DG,Mod,[]);
                   2093:                for ( Ind = [], I = length(G)-1; I >= 0; I-- ) Ind = cons(I,Ind);
                   2094:        }
                   2095:        DF = dp_ptod(F,V); DFI = dp_ptod(1,V);
                   2096:        if ( Mod ) {
                   2097:                DF = dp_mod(DF,Mod,[]); DFI = dp_mod(DFI,Mod,[]);
                   2098:                setmod(Mod);
                   2099:        }
                   2100:        for ( I = 0; I < 3; I++ ) {
                   2101:                DFI = Mod?dp_nf_mod(Ind,DF*DFI,DG,0,Mod):dp_nf(Ind,DF*DFI,DG,0);
                   2102:                if ( !DFI ) return 0;
                   2103:        }
                   2104:        NV = ttttt;
                   2105:        if ( IsGB )
                   2106:                T = nd_gr(append(G,[NV*F-1]),cons(NV,V),Mod,[[0,1],[0,length(V)]]
                   2107:                        |gbblock=[[0,length(G)]]);
                   2108:        else
                   2109:                T = nd_gr(append(G,[NV*F-1]),cons(NV,V),Mod,[[0,1],[0,length(V)]]);
                   2110:        if ( type(car(T)) == 1 ) return 0;
                   2111:        G0 = elimination(T,V);
                   2112:        return G0;
                   2113: }
                   2114:
1.1       noro     2115: def modular_radical_membership(F,G,V) {
                   2116:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2117:        if ( Mod )
                   2118:                return radical_membership(F,G,V|mod=Mod);
                   2119:
                   2120:        F = gen_nf(F,G,V,0,0);
                   2121:        if ( !F ) return 0;
                   2122:        NV = ttttt;
                   2123:        for ( J = 0; ; J++ ) {
                   2124:                Mod = lprime(J);
                   2125:                H = map(dp_hc,map(dp_ptod,G,V));
                   2126:                for ( ; H != []; H = cdr(H) ) if ( !(car(H)%Mod) ) break;
                   2127:                if ( H != [] ) continue;
                   2128:
                   2129:                T = nd_f4(cons(NV*F-1,G),cons(NV,V),Mod,0);
                   2130:                if ( type(car(T)) == 1 ) {
                   2131:                        I = radical_membership_rep(F,G,V,-1,0,Mod);
                   2132:                        I1 = radical_membership_rep(F,G,V,I,0,0);
                   2133:                        if ( I1 > 0 ) return 0;
                   2134:                }
                   2135:                return radical_membership(F,G,V);
                   2136:        }
                   2137: }
                   2138:
                   2139: def radical_membership_rep(F,G,V,Max,Ord,Mod) {
                   2140:        Ft = F;
                   2141:        I = 1;
                   2142:        while ( Max < 0 || I <= Max ) {
                   2143:                Ft = gen_nf(Ft,G,V,Ord,Mod);
                   2144:                if ( !Ft ) return I;
                   2145:                Ft *= F;
                   2146:                I++;
                   2147:        }
                   2148:        return -1;
                   2149: }
                   2150:
                   2151: def ideal_product(A,B,V)
                   2152: {
                   2153:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2154:        dp_ord(0);
                   2155:        DA = map(dp_ptod,A,V);
                   2156:        DB = map(dp_ptod,B,V);
                   2157:        DegA = map(dp_td,DA);
                   2158:        DegB = map(dp_td,DB);
                   2159:        for ( PA = [], T = A, DT = DegA; T != []; T = cdr(T), DT = cdr(DT) )
                   2160:                PA = cons([car(T),car(DT)],PA);
                   2161:        PA = reverse(PA);
                   2162:        for ( PB = [], T = B, DT = DegB; T != []; T = cdr(T), DT = cdr(DT) )
                   2163:                PB = cons([car(T),car(DT)],PB);
                   2164:        PB = reverse(PB);
                   2165:        R = [];
                   2166:        for ( T = PA; T != []; T = cdr(T) )
                   2167:                for ( S = PB; S != []; S = cdr(S) )
                   2168:                        R = cons([car(T)[0]*car(S)[0],car(T)[1]+car(S)[1]],R);
1.9       ohara    2169:        T = qsort(R,noro_pd.comp_by_second);
1.1       noro     2170:        T = map(first,T);
                   2171:        Len = length(A)>length(B)?length(A):length(B);
                   2172:        Len *= 2;
                   2173:        L = sep_list(T,Len); B0 = L[0]; B1 = L[1];
                   2174:        R = fast_gb(B0,V,Mod,0);
                   2175:        while ( B1 != [] ) {
                   2176:                print(length(B1));
                   2177:                L = sep_list(B1,Len);
                   2178:                B0 = L[0]; B1 = L[1];
                   2179:                R = fast_gb(append(R,B0),V,Mod,0|gbblock=[[0,length(R)]],nora=1);
                   2180:        }
                   2181:        return R;
                   2182: }
                   2183:
                   2184: def saturation(GNV,F,V)
                   2185: {
                   2186:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2187:        G = GNV[0]; NV = GNV[1];
                   2188:        if ( Mod )
                   2189:                G1 = nd_gr(cons(NV*F-1,G),cons(NV,V),Mod,[[0,1],[0,length(V)]]);
                   2190:        else if ( Procs ) {
                   2191:                Arg0 = ["nd_gr_trace",
                   2192:                cons(NV*F-1,G),cons(NV,V),0,GBCheck,[[0,1],[0,length(V)]]];
                   2193:                Arg1 = ["nd_gr_trace",
                   2194:                cons(NV*F-1,G),cons(NV,V),1,GBCheck,[[0,1],[0,length(V)]]];
                   2195:                G1 = competitive_exec(Procs,Arg0,Arg1);
                   2196:        } else
                   2197:                G1 = nd_gr(cons(NV*F-1,G),cons(NV,V),0,[[0,1],[0,length(V)]]);
                   2198:        return elimination(G1,V);
                   2199: }
                   2200:
                   2201: def sat(G,F,V)
                   2202: {
                   2203:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2204:        if ( type(IsGB=getopt(isgb)) == -1 ) IsGB = 0;
                   2205:        NV = ttttt;
                   2206:        if ( Mod )
                   2207:                G1 = nd_gr(cons(NV*F-1,G),cons(NV,V),Mod,[[0,1],[0,length(V)]]);
                   2208:        else if ( Procs ) {
                   2209:                Arg0 = ["nd_gr_trace",
                   2210:                cons(NV*F-1,G),cons(NV,V),0,GBCheck,[[0,1],[0,length(V)]]];
                   2211:                Arg1 = ["nd_gr_trace",
                   2212:                cons(NV*F-1,G),cons(NV,V),1,GBCheck,[[0,1],[0,length(V)]]];
                   2213:                G1 = competitive_exec(Procs,Arg0,Arg1);
                   2214:        } else {
                   2215:                B1 = append(G,[NV*F-1]);
                   2216:                V1 = cons(NV,V);
                   2217:                Ord1 = [[0,1],[0,length(V)]];
                   2218:                if ( IsGB )
                   2219:                        G1 = nd_gr(B1,V1,0,Ord1|gbblock=[[0,length(G)]]);
                   2220:                else
                   2221:                        G1 = nd_gr(B1,V1,0,Ord1);
                   2222:        }
                   2223:        return elimination(G1,V);
                   2224: }
                   2225:
1.4       noro     2226: def isat(B,S,V)
                   2227: {
                   2228:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2229:        if ( type(IsGB=getopt(isgb)) == -1 ) IsGB = 0;
                   2230:        F = cdr(fctr(S));
                   2231:        R = B;
                   2232:        for ( T = F; T != []; T = cdr(T) )
                   2233:                R = sat(R,car(T)[0],V|mod=Mod,isgb=IsGB);
                   2234:        return R;
                   2235: }
                   2236:
1.11    ! noro     2237: /* buggy; do not use */
1.1       noro     2238: def satind(G,F,V)
                   2239: {
                   2240:        if ( type(Block=getopt(gbblock)) == -1 ) Block = 0;
                   2241:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2242:        NV = ttttt;
                   2243:        N = length(V);
                   2244:        B = append(G,[NV*F-1]);
                   2245:        V1 = cons(NV,V);
                   2246:        Ord1 = [[0,1],[0,N]];
                   2247:        if ( Mod )
                   2248:                if ( Block )
                   2249:                        D = nd_gr(B,V1,Mod,Ord1|nora=1,gentrace=1,gbblock=Block);
                   2250:                else
                   2251:                        D = nd_gr(B,V1,Mod,Ord1|nora=1,gentrace=1);
                   2252:        else
                   2253:                if ( Block )
                   2254:                        D = nd_gr_trace(B,V1,SatHomo,GBCheck,Ord1
                   2255:                                |nora=1,gentrace=1,gbblock=Block);
                   2256:                else
                   2257:                        D = nd_gr_trace(B,V1,SatHomo,GBCheck,Ord1
                   2258:                                |nora=1,gentrace=1);
                   2259:        G1 = D[0];
                   2260:        Len = length(G1);
                   2261:        Deg = compute_deg(B,V1,NV,D);
                   2262:        D1 = 0;
                   2263:        R = [];
                   2264:        M = length(B);
                   2265:        for ( I = 0; I < Len; I++ ) {
                   2266:                if ( !member(NV,vars(G1[I])) ) {
                   2267:                        for ( J = 1; J < M; J++ )
                   2268:                                D1 = MAX(D1,Deg[I][J]);
                   2269:                        R = cons(G1[I],R);
                   2270:                }
                   2271:        }
                   2272:        return [reverse(R),D1];
                   2273: }
                   2274:
1.11    ! noro     2275: def sat_ind_var(G,F,V)
        !          2276: {
        !          2277:        if ( type(Ord=getopt(ord)) == -1 ) Ord = 0;
        !          2278:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
        !          2279:   V0 = append(setminus(V,[F]),[F]);
        !          2280:   G0 = nd_gr(G,V0,Mod,0);
        !          2281:   M = 0;
        !          2282:   for ( G1 = [], T = G0; T != []; T = cdr(T) ) {
        !          2283:     S = car(T);
        !          2284:     M1 = mindeg(S,F);
        !          2285:     S = sdiv(S,F^M1);
        !          2286:     G1 = cons(S,G1);
        !          2287:     if ( M1 > M ) M = M1;
        !          2288:   }
        !          2289:   G1 = nd_gr(G1,V,Mod,Ord);
        !          2290:   return [G1,M];
        !          2291: }
        !          2292:
1.1       noro     2293: def sat_ind(G,F,V)
                   2294: {
                   2295:        if ( type(Ord=getopt(ord)) == -1 ) Ord = 0;
                   2296:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
1.11    ! noro     2297:   if ( F == var(F) )
        !          2298:     return sat_ind_var(G,F,V|ord=Ord,mod=Mod);
1.1       noro     2299:        NV = ttttt;
                   2300:        F = gen_nf(F,G,V,Ord,Mod);
                   2301:        for ( I = 0, GI = G; ; I++ ) {
                   2302:                G1 = colon(GI,F,V|mod=Mod,ord=Ord);
                   2303:                if ( ideal_inclusion(G1,GI,V,Ord|mod=Mod) )  {
                   2304:                        return [GI,I];
                   2305:                }
                   2306:                else GI = G1;
                   2307:        }
                   2308: }
                   2309:
                   2310: def colon(G,F,V)
                   2311: {
                   2312:        if ( type(Ord=getopt(ord)) == -1 ) Ord = 0;
                   2313:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2314:        if ( type(IsGB=getopt(isgb)) == -1 ) IsGB = 0;
                   2315:        F = gen_nf(F,G,V,Ord,Mod);
                   2316:        if ( !F ) return [1];
                   2317:        if ( IsGB )
                   2318:                T = ideal_intersection(G,[F],V,Ord|gbblock=[[0,length(G)]],mod=Mod);
                   2319:        else
                   2320:                T = ideal_intersection(G,[F],V,Ord|mod=Mod);
1.4       noro     2321:        Gen = Mod?map(sdivm,T,F,Mod):map(ptozp,map(sdiv,T,F));
                   2322:        return nd_gr(Gen,V,Mod,Ord);
1.1       noro     2323: }
                   2324:
                   2325: #if 1
                   2326: def ideal_colon(G,F,V)
                   2327: {
                   2328:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2329:        G = nd_gr(G,V,Mod,0);
                   2330:        C = [1];
                   2331:        TV = ttttt;
1.9       ohara    2332:        F = qsort(F,noro_pd.comp_tdeg);
1.1       noro     2333:        for ( T = F; T != []; T = cdr(T) ) {
                   2334:                S = colon(G,car(T),V|isgb=1,mod=Mod);
                   2335:                if ( type(S[0])!= 1 ) {
                   2336:                        C = nd_gr(append(vtol(ltov(C)*TV),vtol(ltov(S)*(1-TV))),
                   2337:                                cons(TV,V),Mod,[[0,1],[Ord,length(V)]]|gbblock=[[0,length(C)]]);
                   2338:                        C = elimination(C,V);
                   2339:                }
                   2340:        }
                   2341:        return C;
                   2342: }
                   2343: #else
                   2344: def ideal_colon(G,F,V)
                   2345: {
                   2346:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2347:        G = nd_gr(G,V,Mod,0);
                   2348:        for ( T = F, L = []; T != []; T = cdr(T) ) {
                   2349:                C = colon(G,car(T),V|isgb=1,mod=Mod);
                   2350:                if ( type(C[0]) != 1 ) L = cons(C,L);
                   2351:        }
                   2352:        L = reverse(L);
                   2353:        return ideal_list_intersection(L,V,0|mod=Mod);
                   2354: }
                   2355:
                   2356: #endif
                   2357:
                   2358: def member(A,L)
                   2359: {
                   2360:        for ( ; L != []; L = cdr(L) )
                   2361:                if ( car(L) == A ) return 1;
                   2362:        return 0;
                   2363: }
                   2364:
                   2365: def mingen(B,V) {
                   2366:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2367:        Data = nd_gr(B,V,Mod,O|gentrace=1,gensyz=1);
1.4       noro     2368:     G = Data[0]; STrace = Data[6];
1.1       noro     2369:        N = length(G);
1.4       noro     2370:        S = compute_gbsyz(N,V,STrace,Mod);
                   2371:        for ( T = S, R = []; T != []; T = cdr(T) ) {
                   2372:                for ( A = car(T); A1 = dp_rest(A); A = A1);
                   2373:                if ( type(dp_hc(A)) ==1 ) R = cons(dp_etov(A)[0],R);
                   2374:        }
                   2375:        for ( I = 0, U = []; I < N; I++ ) if ( !member(I,R) ) U = cons(G[I],U);
1.1       noro     2376:        return U;
                   2377: }
                   2378:
1.4       noro     2379: def compute_gbsyz(N,V,Trace,Mod)
1.1       noro     2380: {
                   2381:        P = vector(N);
1.4       noro     2382:        for ( I = 0; I < N; I++ ) P[I] = dp_ptod(x^I,[x]);
                   2383:        for ( U = [], T = Trace; T != []; T = cdr(T) ) {
1.1       noro     2384:                Ti = car(T);
                   2385:                if ( Ti[0] != -1 ) error("Input is not a GB");
1.4       noro     2386:                R = recompute_trace(Ti[1],P,V,Mod);
                   2387:                U = cons(R,U);
1.1       noro     2388:        }
                   2389:        return reverse(U);
                   2390: }
                   2391:
1.4       noro     2392: def recompute_trace(Ti,P,V,Mod)
1.1       noro     2393: {
                   2394:   for ( Num = 0, Den = 1; Ti != []; Ti = cdr(Ti) ) {
1.4       noro     2395:     Sj = car(Ti); Dj = Sj[0]; Ij =Sj[1]; Mj = dp_dtop(Sj[2],V); Cj = Sj[3];
1.1       noro     2396:     /* Num/Den <- (Dj*Num+Den*Mj*P[Ij])/(Den*Cj) */
1.4       noro     2397:     if ( Dj ) Num = (Dj*Num+Den*Mj*P[Ij]);
1.1       noro     2398:     Den *= Cj;
                   2399:   }
1.4       noro     2400:   return Num;
1.1       noro     2401: }
                   2402:
                   2403: def ideal_sat(G,F,V)
                   2404: {
                   2405:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2406:        G = nd_gr(G,V,Mod,0);
                   2407:        for ( T = F, L = []; T != []; T = cdr(T) )
                   2408:                L = cons(sat(G,car(T),V|mod=Mod),L);
                   2409:        L = reverse(L);
                   2410:        return ideal_list_intersection(L,V,0|mod=Mod);
                   2411: }
                   2412:
                   2413: def ideal_inclusion(F,G,V,O)
                   2414: {
                   2415:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2416:        for ( T = F; T != []; T = cdr(T) )
                   2417:                if ( gen_nf(car(T),G,V,O,Mod) ) return 0;
                   2418:        return 1;
                   2419: }
                   2420:
                   2421: /* remove redundant components */
                   2422:
                   2423: def qd_simp_comp(QP,V)
                   2424: {
                   2425:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2426:        R = ltov(QP);
                   2427:        N = length(R);
                   2428:        for ( I = 0; I < N; I++ ) {
                   2429:                if ( R[I] ) {
                   2430:                        QI = R[I][0]; PI = R[I][1];
                   2431:                        for ( J = I+1; J < N; J++ )
                   2432:                                if ( R[J] && gen_gb_comp(PI,R[J][1],Mod) ) {
                   2433:                                        QI = ideal_intersection(QI,R[J][0],V,0|mod=Mod);
                   2434:                                        R[J] = 0;
                   2435:                                }
                   2436:                        R[I] = [QI,PI];
                   2437:                }
                   2438:        }
                   2439:        for ( I = N-1, S = []; I >= 0; I-- )
                   2440:                if ( R[I] ) S = cons(R[I],S);
                   2441:        return S;
                   2442: }
                   2443:
                   2444: def qd_remove_redundant_comp(G,Iso,Emb,V,Ord)
                   2445: {
                   2446:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2447:        IsoInt = ideal_list_intersection(map(first,Iso),V,Ord|mod=Mod);
                   2448:        Emb = qd_simp_comp(Emb,V|mod=Mod);
                   2449:        Emb = reverse(qsort(Emb));
                   2450:        A = ltov(Emb); N = length(A);
                   2451:        Pre = IsoInt; Post = vector(N+1);
                   2452:        for ( Post[N] = IsoInt, I = N-1; I >= 1; I-- )
                   2453:                Post[I] = ideal_intersection(Post[I+1],A[I][0],V,Ord|mod=Mod);
                   2454:        for ( I = 0; I < N; I++ ) {
                   2455:                print(".",2);
                   2456:                Int = ideal_intersection(Pre,Post[I+1],V,Ord|mod=Mod);
                   2457:                if ( gen_gb_comp(Int,G,Mod) ) A[I] = 0;
                   2458:                else
                   2459:                        Pre = ideal_intersection(Pre,A[I][0],V,Ord|mod=Mod);
                   2460:        }
                   2461:        for ( T = [], I = 0; I < N; I++ )
                   2462:                if ( A[I] ) T = cons(A[I],T);
                   2463:        return reverse(T);
                   2464: }
                   2465:
                   2466: def pd_simp_comp(PL,V)
                   2467: {
                   2468:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2469:        if ( type(First=getopt(first)) == -1 ) First = 0;
                   2470:        A = ltov(PL); N = length(A);
                   2471:        if ( N == 1 )  return PL;
                   2472:        for ( I = 0; I < N; I++ ) {
                   2473:                if ( !A[I] ) continue;
                   2474:                AI = First?A[I][0]:A[I];
                   2475:                for ( J = 0; J < N; J++ ) {
                   2476:                        if ( J == I || !A[J] ) continue;
                   2477:                        AJ = First?A[J][0]:A[J];
                   2478:                        if ( gen_gb_comp(AI,AJ,Mod) || ideal_inclusion(AI,AJ,V,Ord|mod=Mod) )
                   2479:                                A[J] = 0;
                   2480:                }
                   2481:        }
                   2482:        for ( I = 0, T = []; I < N; I++ ) if ( A[I] ) T = cons(A[I],T);
                   2483:        return reverse(T);
                   2484: }
                   2485:
                   2486: def pd_remove_redundant_comp(G,P,V,Ord)
                   2487: {
                   2488:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2489:        if ( type(First=getopt(first)) == -1 ) First = 0;
                   2490:        if ( length(P) == 1 )  return P;
                   2491:
                   2492:        A = ltov(P); N = length(A);
                   2493:        for ( I = 0; I < N; I++ ) {
                   2494:                if ( !A[I] ) continue;
                   2495:                for ( J = I+1; J < N; J++ )
                   2496:                        if ( A[J] &&
                   2497:                                gen_gb_comp(First?A[I][0]:A[I],First?A[J][0]:A[J],Mod) ) A[J] = 0;
                   2498:        }
                   2499:        for ( I = 0, T = []; I < N; I++ ) if ( A[I] ) T = cons(A[I],T);
                   2500:        A = ltov(reverse(T)); N = length(A);
                   2501:        Pre = [1]; Post = vector(N+1);
                   2502:        for ( Post[N] = [1], I = N-1; I >= 1; I-- )
                   2503:                Post[I] = ideal_intersection(Post[I+1],First?A[I][0]:A[I],V,Ord|mod=Mod);
                   2504:        for ( I = 0; I < N; I++ ) {
                   2505:                Int = ideal_intersection(Pre,Post[I+1],V,Ord|mod=Mod);
                   2506:                if ( gen_gb_comp(Int,G,Mod) ) A[I] = 0;
                   2507:                else
                   2508:                        Pre = ideal_intersection(Pre,First?A[I][0]:A[I],V,Ord|mod=Mod);
                   2509:        }
                   2510:        for ( T = [], I = 0; I < N; I++ ) if ( A[I] ) T = cons(A[I],T);
                   2511:        return reverse(T);
                   2512: }
                   2513:
1.7       noro     2514: def remove_identical_comp(L)
                   2515: {
                   2516:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2517:        if ( length(L) == 1 )  return L;
                   2518:
                   2519:        A = ltov(L); N = length(A);
                   2520:        for ( I = 0; I < N; I++ ) {
                   2521:                if ( !A[I] ) continue;
                   2522:                for ( J = I+1; J < N; J++ )
                   2523:                        if ( A[J] &&
                   2524:                                gen_gb_comp(A[I],A[J],Mod) ) A[J] = 0;
                   2525:        }
                   2526:        for ( I = 0, T = []; I < N; I++ ) if ( A[I] ) T = cons(A[I],T);
                   2527:        return reverse(T);
                   2528: }
                   2529:
1.1       noro     2530: /* polynomial operations */
                   2531:
                   2532: def ppart(F,V,Mod)
                   2533: {
                   2534:        if ( !Mod )
                   2535:                G = nd_gr([F],[V],0,0);
                   2536:        else
                   2537:                G = dp_gr_mod_main([F],[V],0,Mod,0);
                   2538:        return G[0];
                   2539: }
                   2540:
                   2541:
                   2542: def sq(F,Mod)
                   2543: {
                   2544:        if ( !F ) return 0;
                   2545:        A = cdr(gen_fctr(F,Mod));
                   2546:        for ( R = 1; A != []; A = cdr(A) )
                   2547:                R *= car(car(A));
                   2548:        return R;
                   2549: }
                   2550:
                   2551: def lcfactor(G,V,O,Mod)
                   2552: {
                   2553:        O0 = dp_ord(); dp_ord(O);
                   2554:        C = [];
                   2555:        for ( T = G; T != []; T = cdr(T) ) {
                   2556:                C1 = dp_hc(dp_ptod(car(T),V));
                   2557:                S = gen_fctr(C1,Mod);
                   2558:                for ( S = cdr(S); S != []; S = cdr(S) )
                   2559:                        if ( !member(S[0][0],C) ) C = cons(S[0][0],C);
                   2560:        }
                   2561:        dp_ord(O0);
                   2562:        return C;
                   2563: }
                   2564:
                   2565: def gen_fctr(F,Mod)
                   2566: {
                   2567:        if ( Mod ) return modfctr(F,Mod);
                   2568:        else return fctr(F);
                   2569: }
                   2570:
                   2571: def gen_mptop(F)
                   2572: {
                   2573:        if ( !F ) return F;
                   2574:        else if ( type(F)==1 )
                   2575:                if ( ntype(F)==5 ) return mptop(F);
                   2576:                else return F;
                   2577:        else {
                   2578:                V = var(F);
                   2579:                D = deg(F,V);
                   2580:                for ( R = 0, I = 0; I <= D; I++ )
                   2581:                        if ( C = coef(F,I,V) ) R += gen_mptop(C)*V^I;
                   2582:                return R;
                   2583:        }
                   2584: }
                   2585:
                   2586: def gen_nf(F,G,V,Ord,Mod)
                   2587: {
                   2588:        if ( !Mod ) return p_nf(F,G,V,Ord);
                   2589:
                   2590:        setmod(Mod);
                   2591:        dp_ord(Ord); DF = dp_mod(dp_ptod(F,V),Mod,[]);
                   2592:        N = length(G); DG = newvect(N);
                   2593:        for ( I = N-1, IL = []; I >= 0; I-- ) {
                   2594:                DG[I] = dp_mod(dp_ptod(G[I],V),Mod,[]);
                   2595:                IL = cons(I,IL);
                   2596:        }
                   2597:        T = dp_nf_mod(IL,DF,DG,1,Mod);
                   2598:        for ( R = 0; T; T = dp_rest(T) )
                   2599:                R += gen_mptop(dp_hc(T))*dp_dtop(dp_ht(T),V);
                   2600:        return R;
                   2601: }
                   2602:
                   2603: /* Ti = [D,I,M,C] */
                   2604:
                   2605: def compute_deg0(Ti,P,V,TV)
                   2606: {
                   2607:        N = length(P[0]);
                   2608:        Num = vector(N);
                   2609:        for ( I = 0; I < N; I++ ) Num[I] = -1;
                   2610:        for ( ; Ti != []; Ti = cdr(Ti) ) {
                   2611:                Sj = car(Ti);
                   2612:                Dj = Sj[0];
                   2613:                Ij =Sj[1];
                   2614:                Mj = deg(type(Sj[2])==9?dp_dtop(Sj[2],V):Sj[2],TV);
                   2615:                Pj = P[Ij];
                   2616:                if ( Dj )
                   2617:                        for ( I = 0; I < N; I++ )
                   2618:                                if ( Pj[I] >= 0 ) {
                   2619:                                        T = Mj+Pj[I];
                   2620:                                        Num[I] = MAX(Num[I],T);
                   2621:                                }
                   2622:        }
                   2623:        return Num;
                   2624: }
                   2625:
                   2626: def compute_deg(B,V,TV,Data)
                   2627: {
                   2628:        GB = Data[0];
                   2629:        Homo = Data[1];
                   2630:        Trace = Data[2];
                   2631:        IntRed = Data[3];
                   2632:        Ind = Data[4];
                   2633:        DB = map(dp_ptod,B,V);
                   2634:        if ( Homo ) {
                   2635:                DB = map(dp_homo,DB);
                   2636:                V0 = append(V,[hhh]);
                   2637:        } else
                   2638:                V0 = V;
                   2639:        Perm = Trace[0]; Trace = cdr(Trace);
                   2640:        for ( I = length(Perm)-1, T = Trace; T != []; T = cdr(T) )
                   2641:                if ( (J=car(T)[0]) > I ) I = J;
                   2642:        N = I+1;
                   2643:        N0 = length(B);
                   2644:        P = vector(N);
                   2645:        for ( T = Perm, I = 0; T != []; T = cdr(T), I++ ) {
                   2646:                Pi = car(T);
                   2647:                C = vector(N0);
                   2648:                for ( J = 0; J < N0; J++ ) C[J] = -1;
                   2649:                C[Pi[1]] = 0;
                   2650:                P[Pi[0]] = C;
                   2651:        }
                   2652:        for ( T = Trace; T != []; T = cdr(T) ) {
                   2653:                Ti = car(T); P[Ti[0]] = compute_deg0(Ti[1],P,V0,TV);
                   2654:        }
                   2655:        M = length(Ind);
                   2656:        for ( T = IntRed; T != []; T = cdr(T) ) {
                   2657:                Ti = car(T); P[Ti[0]] = compute_deg0(Ti[1],P,V,TV);
                   2658:        }
                   2659:        R = [];
                   2660:        for ( J = 0; J < M; J++ ) {
                   2661:                U = P[Ind[J]];
                   2662:                R = cons(U,R);
                   2663:        }
                   2664:        return reverse(R);
                   2665: }
                   2666:
                   2667: /* set theoretic functions */
                   2668:
                   2669: def member(A,S)
                   2670: {
                   2671:        for ( ; S != []; S = cdr(S) )
                   2672:                if ( car(S) == A ) return 1;
                   2673:        return 0;
                   2674: }
                   2675:
                   2676: def elimination(G,V) {
                   2677:        for ( R = [], T = G; T != []; T = cdr(T) )
                   2678:                if ( setminus(vars(car(T)),V) == [] ) R =cons(car(T),R);
                   2679:        return R;
                   2680: }
                   2681:
                   2682: def setintersection(A,B)
                   2683: {
                   2684:        for ( L = []; A != []; A = cdr(A) )
                   2685:                if ( member(car(A),B) )
                   2686:                        L = cons(car(A),L);
                   2687:        return L;
                   2688: }
                   2689:
                   2690: def setminus(A,B) {
                   2691:        for ( T = reverse(A), R = []; T != []; T = cdr(T) ) {
                   2692:                for ( S = B, M = car(T); S != []; S = cdr(S) )
                   2693:                        if ( car(S) == M ) break;
                   2694:                if ( S == [] ) R = cons(M,R);
                   2695:        }
                   2696:        return R;
                   2697: }
                   2698:
                   2699: def sep_list(L,N)
                   2700: {
                   2701:        if ( length(L) <= N ) return [L,[]];
                   2702:        R = [];
                   2703:        for ( T = L, I = 0; I < N; I++, T = cdr(T) )
                   2704:                R = cons(car(T),R);
                   2705:        return [reverse(R),T];
                   2706: }
                   2707:
                   2708: def first(L)
                   2709: {
                   2710:        return L[0];
                   2711: }
                   2712:
                   2713: def second(L)
                   2714: {
                   2715:        return L[1];
                   2716: }
                   2717:
                   2718: def third(L)
                   2719: {
                   2720:        return L[2];
                   2721: }
                   2722:
                   2723: def first_second(L)
                   2724: {
                   2725:        return [L[0],L[1]];
                   2726: }
                   2727:
                   2728: def comp_tord(A,B)
                   2729: {
                   2730:        DA = dp_ht(A);
                   2731:        DB = dp_ht(B);
                   2732:        if ( DA > DB ) return 1;
                   2733:        else if ( DA < DB ) return -1;
                   2734:        else return 0;
                   2735: }
                   2736:
                   2737: def comp_tdeg(A,B)
                   2738: {
                   2739:        DA = tdeg(A);
                   2740:        DB = tdeg(B);
                   2741:        if ( DA > DB ) return 1;
                   2742:        else if ( DA < DB ) return -1;
                   2743:        else return 0;
                   2744: }
                   2745:
                   2746: def comp_tdeg_first(A,B)
                   2747: {
                   2748:        DA = tdeg(A[0]);
                   2749:        DB = tdeg(B[0]);
                   2750:        if ( DA > DB ) return 1;
                   2751:        else if ( DA < DB ) return -1;
                   2752:        else return 0;
                   2753: }
                   2754:
                   2755: def comp_third_tdeg(A,B)
                   2756: {
                   2757:        if ( A[2] > B[2] ) return 1;
                   2758:        if ( A[2] < B[2] ) return -1;
                   2759:        DA = tdeg(A[0]);
                   2760:        DB = tdeg(B[0]);
                   2761:        if ( DA > DB ) return 1;
                   2762:        else if ( DA < DB ) return -1;
                   2763:        else return 0;
                   2764: }
                   2765:
                   2766: def tdeg(P)
                   2767: {
                   2768:        dp_ord(0);
                   2769:        return dp_td(dp_ptod(P,vars(P)));
                   2770: }
                   2771:
                   2772: def comp_by_ord(A,B)
                   2773: {
                   2774:        if ( dp_ht(A) > dp_ht(B) ) return 1;
                   2775:        else if ( dp_ht(A) < dp_ht(B) ) return -1;
                   2776:        else return 0;
                   2777: }
                   2778:
                   2779: def comp_by_second(A,B)
                   2780: {
                   2781:        if ( A[1] > B[1] ) return 1;
                   2782:        else if ( A[1] < B[1] ) return -1;
                   2783:        else return 0;
                   2784: }
                   2785:
                   2786: def get_lc(F)
                   2787: {
                   2788:        if ( type(F)==1 ) return F;
                   2789:        V = var(F);
                   2790:        D = deg(F,V);
                   2791:        return get_lc(coef(F,D,V));
                   2792: }
                   2793:
                   2794: def tomonic(F,Mod)
                   2795: {
                   2796:        C = get_lc(F);
                   2797:        IC = inv(C,Mod);
                   2798:        return (IC*F)%Mod;
                   2799: }
                   2800:
                   2801: def gen_gb_comp(A,B,Mod)
                   2802: {
                   2803:        if ( !Mod ) return gb_comp(A,B);
                   2804:        LA = length(A); LB = length(B);
                   2805:        if ( LA != LB ) return 0;
                   2806:        A = map(tomonic,A,Mod);
                   2807:        B = map(tomonic,B,Mod);
                   2808:        A = qsort(A); B = qsort(B);
                   2809:        if ( A != B ) return 0;
                   2810:        return 1;
                   2811: }
                   2812:
                   2813: def prod(L)
                   2814: {
                   2815:        for ( R = 1; L != []; L = cdr(L) )
                   2816:                R *= car(L);
                   2817:        return R;
                   2818: }
                   2819:
                   2820: def monodec0(B,V)
                   2821: {
                   2822:        M = monodec(B,V);
                   2823:        return map(vars,M);
                   2824: }
                   2825:
                   2826: def monodec(B,V)
                   2827: {
                   2828:        B = map(sq,B,0);
                   2829:        G = nd_gr_postproc(B,V,0,0,0);
                   2830:        V = vars(G);
                   2831:        N = length(V);
                   2832:        if ( N == 0 ) return G == [] ? [[]] : [];
                   2833:        if ( N == 1 ) return G;
                   2834:        if ( N < 20 ) {
                   2835:                T = dp_mono_raddec(G,V);
                   2836:                return map(prod,T);
                   2837:        }
                   2838:        X = car(V); W = cdr(V);
                   2839:        D0 = monodec(map(subst,B,X,0),W);
                   2840:        T0 = map(dp_ptod,D0,W);
                   2841:        D1 = monodec(map(subst,B,X,1),W);
                   2842:        T1 = map(dp_ptod,D1,W);
1.4       noro     2843: #if 0
1.1       noro     2844:        for ( T = T1; T != []; T = cdr(T) ) {
                   2845:                for ( M = car(T), S1 = [], S = T0; S != []; S = cdr(S) )
                   2846:                        if ( !dp_redble(car(S),M) ) S1= cons(car(S),S1);
                   2847:                T0 = S1;
                   2848:        }
1.4       noro     2849: #else
                   2850:        T0 = dp_mono_reduce(T0,T1);
                   2851: #endif
1.1       noro     2852:        D0 = map(dp_dtop,T0,W);
                   2853:        D0 = vtol(X*ltov(D0));
                   2854:        return append(D0,D1);
                   2855: }
                   2856:
                   2857: def separator(P,V)
                   2858: {
                   2859:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
                   2860:        N = length(P);
1.4       noro     2861:        dp_ord(0);
                   2862:        DP = vector(N);
1.9       ohara    2863:        for ( I = 0; I < N; I++ ) DP[I] = qsort(ltov(map(dp_ptod,P[I][0],V)),noro_pd.comp_tord);
1.4       noro     2864:        if ( Mod )
                   2865:                for ( I = 0; I < N; I++ ) DP[I] = map(dp_mod,DP[I],Mod,[]);
                   2866:        Ind = vector(N);
                   2867:        for ( I = 0; I < N; I++ ) {
                   2868:                for ( K = [], J = length(DP[I])-1; J >= 0; J-- ) K = cons(J,K);
                   2869:                Ind[I] = K;
                   2870:        }
                   2871:        S = vector(N);
                   2872:        for ( I = 0; I < N; I++ ) S[I] = 1;
1.1       noro     2873:        for ( I = 0; I < N; I++ ) {
1.4       noro     2874:                print(".",2);
1.1       noro     2875:                for ( J = 0; J < N; J++ ) {
                   2876:                        if ( J == I ) continue;
1.4       noro     2877:                        T = DP[I]; L = length(T);
                   2878:                        if ( Mod ) {
                   2879:                                for ( K = 0; K < L; K++ )
                   2880:                                        if ( dp_nf_mod(Ind[J],T[K],DP[J],0,Mod) ) break;
                   2881:                        } else {
                   2882:                                for ( K = 0; K < L; K++ )
                   2883:                                        if ( dp_nf(Ind[J],T[K],DP[J],0) ) break;
                   2884:                        }
                   2885:                        S[J] = lcm(S[J],dp_dtop(T[K],V));
1.1       noro     2886:                }
                   2887:        }
1.4       noro     2888:        print("");
1.1       noro     2889:        return S;
                   2890: }
                   2891:
                   2892: def prepost(PL,V)
1.4       noro     2893: {
                   2894:        if ( type(Mod=getopt(mod)) == -1 ) Mod = 0;
1.1       noro     2895:        A = ltov(PL); N = length(A);
                   2896:        Pre = vector(N);
                   2897:        Post = vector(N);
                   2898:        R = vector(N);
                   2899:        Pre[0] = [1];
                   2900:        print("pre ",2);
                   2901:        for ( I = 1; I < N; I++, print(".",2) )
1.4       noro     2902:                Pre[I] = ideal_intersection_m(Pre[I-1],A[I-1],V,0|mod=Mod);
1.1       noro     2903:        print("done");
                   2904:        print("post ",2);
                   2905:        Post[N-1] = [1];
                   2906:        for ( I = N-2; I >= 0; I--, print(".",2) )
1.4       noro     2907:                Post[I] = ideal_intersection_m(Post[I+1],A[I+1],V,0|mod=Mod);
1.1       noro     2908:        print("done");
                   2909:        print("int ",2);
                   2910:        for ( I = 0; I < N; I++, print(".",2) )
1.4       noro     2911:                R[I] = ideal_intersection_m(Pre[I],Post[I],V,0|mod=Mod);
1.1       noro     2912:        print("done");
                   2913:        return R;
                   2914: }
                   2915:
                   2916: /* XXX */
                   2917:
                   2918: def call_func(Arg)
                   2919: {
                   2920:        F = car(Arg);
                   2921:        return call(strtov(F),cdr(Arg));
                   2922: }
                   2923:
1.8       noro     2924: def call_func_serial(Arg,Serial)
                   2925: {
                   2926:        F = car(Arg);
                   2927:        return [call(strtov(F),cdr(Arg)),Serial];
                   2928: }
                   2929:
1.1       noro     2930: def competitive_exec(P,Arg0,Arg1)
                   2931: {
                   2932:        P0 = P[0]; P1 = P[1];
                   2933:        ox_cmo_rpc(P0,"noro_pd.call_func",Arg0|sync=1);
                   2934:        ox_cmo_rpc(P1,"noro_pd.call_func",Arg1|sync=1);
                   2935:        F = ox_select(P);
                   2936:        R = ox_get(F[0]);
                   2937:        if ( length(F) == 2 ) {
                   2938:                ox_get(F[1]);
                   2939:        } else {
                   2940:                U = setminus(P,F);
                   2941:                ox_reset(U[0]);
                   2942:        }
                   2943:        return R;
                   2944: }
                   2945:
                   2946:
                   2947: def nd_gr_rat(B,V,PV,Ord1,Ord)
                   2948: {
                   2949:        G = nd_gr(B,append(V,PV),0,Ord1);
                   2950:        G1 = nd_gr_postproc(G,V,0,Ord,0);
                   2951:        return G1;
                   2952: }
                   2953:
                   2954: /* Task[i] = [fname,[arg0,...,argn]] */
                   2955:
                   2956: def para_exec(Proc,Task) {
                   2957:        Free = Proc;
                   2958:        N = length(Task);
                   2959:        R = [];
1.8       noro     2960:        print([N],2); print("->",2);
                   2961:        Serial = 0;
1.1       noro     2962:        while ( N ) {
                   2963:                while ( Task != [] && Free != [] ) {
                   2964:                        T = car(Task); Task = cdr(Task);
1.8       noro     2965:                        ox_rpc(car(Free),"noro_pd.call_func_serial",T,Serial++);
1.1       noro     2966:                        ox_push_cmd(car(Free),258); Free = cdr(Free);
                   2967:                }
                   2968:                Finish0 = Finish = ox_select(Proc);
                   2969:                for ( ; Finish != []; Finish = cdr(Finish) ) {
                   2970:                        print(".",2);
                   2971:                        L = ox_get(car(Finish));
                   2972:                        R = cons(L,R);
                   2973:                        N--;
                   2974:                }
1.8       noro     2975:                print([N],2);
1.1       noro     2976:                Free = append(Free,Finish0);
                   2977:        }
                   2978:        print("");
1.9       ohara    2979:        R = qsort(R,noro_pd.comp_by_second);
1.8       noro     2980:        R = map(first,R);
                   2981:        return R;
1.1       noro     2982: }
1.4       noro     2983:
                   2984: def redbase(B,V,Mod,Ord)
                   2985: {
                   2986:        M = nd_gr_postproc(B,V,Mod,Ord,0);
                   2987:        dp_ord(Ord);
                   2988:        DM = ltov(map(dp_ptod,M,V));
                   2989:        if ( Mod ) DM = map(dp_mod,DM,Mod,[]);
                   2990:        N = length(DM);
                   2991:        for ( Ind = [], I = N-1; I >= 0; I-- ) Ind = cons(I,Ind);
                   2992:        for ( T = B, R = vtol(DM); T != []; T = cdr(T) ) {
                   2993:                D = dp_ptod(car(T),V);
                   2994:                if ( Mod ) D = dp_mod(D,Mod,[]);
                   2995:                D = Mod?dp_nf_mod(Ind,D,DM,1,Mod):dp_nf(Ind,D,DM,1);
                   2996:                if ( D ) R = cons(D,R);
                   2997:        }
1.9       ohara    2998:        D = qsort(R,noro_pd.comp_tord);
1.4       noro     2999:        return map(dp_dtop,D,V);
                   3000: }
                   3001:
                   3002: def witness(A,B,V)
                   3003: {
                   3004:        G = nd_gr(B,V,0,Mod);
                   3005:        L = length(A);
                   3006:        QL = []; PL = [];
                   3007:        for ( I = L-1; I >= 0; I-- ) {
                   3008:                QL = append(map(first,A[I]),QL);
                   3009:                PL = append(map(second,A[I]),PL);
                   3010:        }
                   3011:        N = length(QL);
                   3012:        Qhat = prepost(QL,V);
                   3013:        for ( I = 0, W = []; I < N; I++ ) {
                   3014:                for ( T = Qhat[I]; T != []; T = cdr(T) )
                   3015:                        if ( gen_nf(car(T),QL[I],V,0,Mod) ) break;
                   3016:                Ai = car(T);
                   3017:                Ji = colon(G,Ai,V|isgb=1,mod=Mod);
                   3018:                Ji = nd_gr(Ji,V,Mod,0);
                   3019:                if ( gen_gb_comp(Ji,PL[I],Mod) ) Bi = 1;
                   3020:                else {
                   3021:                        Ki = ideal_colon(Ji,PL[I],V|mod=Mod);
                   3022:                        for ( T = Ki; T != []; T = cdr(T) )
                   3023:                                if ( gen_nf(car(T),Ji,V,0,Mod) ) break;
                   3024:                        Bi = car(T);
                   3025:                }
                   3026:                W = cons(Ai*Bi,W);
                   3027:                Li = colon(G,W[0],V|isgb=1,mod=Mod);
                   3028:                Li = nd_gr(Li,V,Mod,0);
                   3029:                if ( !gen_gb_comp(Li,PL[I],Mod) )
                   3030:                        error("afo");
                   3031:        }
                   3032:        return reverse(W);
                   3033: }
1.1       noro     3034: endmodule$
                   3035: end$

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