=================================================================== RCS file: /home/cvs/OpenXM_contrib2/asir2000/engine/nd.c,v retrieving revision 1.203 retrieving revision 1.212 diff -u -p -r1.203 -r1.212 --- OpenXM_contrib2/asir2000/engine/nd.c 2013/01/31 01:13:47 1.203 +++ OpenXM_contrib2/asir2000/engine/nd.c 2013/09/27 02:35:15 1.212 @@ -1,4 +1,4 @@ -/* $OpenXM: OpenXM_contrib2/asir2000/engine/nd.c,v 1.202 2013/01/30 08:03:18 noro Exp $ */ +/* $OpenXM: OpenXM_contrib2/asir2000/engine/nd.c,v 1.211 2013/09/26 08:55:11 noro Exp $ */ #include "nd.h" @@ -54,6 +54,8 @@ static NODE nd_tracelist; static NODE nd_alltracelist; static int nd_gentrace,nd_gensyz,nd_nora,nd_newelim,nd_intersect; static int *nd_gbblock; +static NODE nd_nzlist,nd_check_splist; +static int nd_splist; NumberField get_numberfield(); UINT *nd_det_compute_bound(NDV **dm,int n,int j); @@ -74,6 +76,7 @@ void conv_ilist(int demand,int trace,NODE g,int **indp void parse_nd_option(NODE opt); void dltondl(int n,DL dl,UINT *r); DP ndvtodp(int mod,NDV p); +DP ndtodp(int mod,ND p); extern int Denominator,DP_Multiple; @@ -1866,6 +1869,29 @@ int do_diagonalize(int sugar,int m) return 1; } +LIST compute_splist() +{ + NODE g,tn0,tn,node; + LIST l0; + ND_pairs d,t; + int i; + Q i1,i2; + + g = 0; d = 0; + for ( i = 0; i < nd_psn; i++ ) { + d = update_pairs(d,g,i,0); + g = update_base(g,i); + } + for ( t = d, tn0 = 0; t; t = NEXT(t) ) { + NEXTNODE(tn0,tn); + STOQ(t->i1,i1); STOQ(t->i2,i2); + node = mknode(2,i1,i2); MKLIST(l0,node); + BDY(tn) = l0; + } + if ( tn0 ) NEXT(tn) = 0; MKLIST(l0,tn0); + return l0; +} + /* return value = 0 => input is not a GB */ NODE nd_gb(int m,int ishomo,int checkonly,int gensyz,int **indp) @@ -1973,6 +1999,45 @@ again: return g; } +/* splist = [[i1,i2],...] */ + +int check_splist(int m,NODE splist) +{ + NODE t,p; + ND_pairs d,r,l; + int stat; + ND h,nf; + + for ( d = 0, t = splist; t; t = NEXT(t) ) { + p = BDY((LIST)BDY(t)); + NEXTND_pairs(d,r); + r->i1 = QTOS((Q)ARG0(p)); r->i2 = QTOS((Q)ARG1(p)); + ndl_lcm(DL(nd_psh[r->i1]),DL(nd_psh[r->i2]),r->lcm); + SG(r) = TD(LCM(r)); /* XXX */ + } + if ( d ) NEXT(r) = 0; + + while ( d ) { +again: + l = nd_minp(d,&d); + stat = nd_sp(m,0,l,&h); + if ( !stat ) { + NEXT(l) = d; d = l; + d = nd_reconstruct(0,d); + goto again; + } + stat = nd_nf(m,0,h,nd_ps,!Top,0,&nf); + if ( !stat ) { + NEXT(l) = d; d = l; + d = nd_reconstruct(0,d); + goto again; + } else if ( nf ) return 0; + if ( DP_Print) { printf("."); fflush(stdout); } + } + if ( DP_Print) { printf("done.\n"); fflush(stdout); } + return 1; +} + int do_diagonalize_trace(int sugar,int m) { int i,nh,stat; @@ -2822,7 +2887,7 @@ void nd_gr(LIST f,LIST v,int m,int homo,int retdp,int Obj obj; NumberField nf; struct order_spec *ord1; - NODE tr,tl1,tl2,tl3,tl4; + NODE tr,tl1,tl2,tl3,tl4,nzlist; LIST l1,l2,l3,l4,l5; int j; Q jq,bpe; @@ -2882,7 +2947,7 @@ void nd_gr(LIST f,LIST v,int m,int homo,int retdp,int max = MAX(e,max); } } - nd_setup_parameters(nvar,max); + nd_setup_parameters(nvar,nd_nzlist?0:max); obpe = nd_bpe; oadv = nmv_adv; oepos = nd_epos; ompos = nd_mpos; ishomo = 1; for ( fd0 = 0, t = BDY(f); t; t = NEXT(t) ) { @@ -2915,10 +2980,19 @@ void nd_gr(LIST f,LIST v,int m,int homo,int retdp,int ndv_homogenize((NDV)BDY(t),obpe,oadv,oepos,ompos); } - ndv_setup(m,0,fd0,nd_gbblock?1:0,0); + ndv_setup(m,0,fd0,(nd_gbblock||nd_splist||nd_check_splist)?1:0,0); if ( nd_gentrace ) { MKLIST(l1,nd_tracelist); MKNODE(nd_alltracelist,l1,0); } + if ( nd_splist ) { + *rp = compute_splist(); + return; + } + if ( nd_check_splist ) { + if ( check_splist(m,nd_check_splist) ) *rp = (LIST)ONE; + else *rp = 0; + return; + } x = f4?nd_f4(m,&perm):nd_gb(m,ishomo || homo,0,0,&perm); if ( !x ) { *rp = 0; return; @@ -2938,11 +3012,12 @@ void nd_gr(LIST f,LIST v,int m,int homo,int retdp,int conv_ilist(nd_demand,0,x,0); goto FINAL; } + if ( nd_gentrace && f4 ) { nzlist = nd_alltracelist; } x = ndv_reducebase(x,perm); - if ( nd_gentrace ) { tl1 = nd_alltracelist; nd_alltracelist = 0; } + if ( nd_gentrace && !f4 ) { tl1 = nd_alltracelist; nd_alltracelist = 0; } x = ndv_reduceall(m,x); cbpe = nd_bpe; - if ( nd_gentrace ) { + if ( nd_gentrace && !f4 ) { tl2 = nd_alltracelist; nd_alltracelist = 0; ndv_check_membership(m,fd0,obpe,oadv,oepos,x); tl3 = nd_alltracelist; nd_alltracelist = 0; @@ -2961,29 +3036,34 @@ FINAL: else BDY(r) = ndvtop(m,CO,vv,BDY(t)); } if ( r0 ) NEXT(r) = 0; - if ( nalg ) + if ( !m && nd_nalg ) r0 = postprocess_algcoef(av,alist,r0); MKLIST(*rp,r0); if ( nd_gentrace ) { - tl1 = reverse_node(tl1); tl2 = reverse_node(tl2); - tl3 = reverse_node(tl3); - /* tl2 = [[i,[[*,j,*,*],...]],...] */ - for ( t = tl2; t; t = NEXT(t) ) { - /* s = [i,[*,j,*,*],...] */ - s = BDY((LIST)BDY(t)); - j = perm[QTOS((Q)ARG0(s))]; STOQ(j,jq); ARG0(s) = (pointer)jq; - for ( s = BDY((LIST)ARG1(s)); s; s = NEXT(s) ) { - j = perm[QTOS((Q)ARG1(BDY((LIST)BDY(s))))]; STOQ(j,jq); - ARG1(BDY((LIST)BDY(s))) = (pointer)jq; + if ( f4 ) { + STOQ(16,bpe); + tr = mknode(4,*rp,(!ishomo&&homo)?ONE:0,BDY(nzlist),bpe); MKLIST(*rp,tr); + } else { + tl1 = reverse_node(tl1); tl2 = reverse_node(tl2); + tl3 = reverse_node(tl3); + /* tl2 = [[i,[[*,j,*,*],...]],...] */ + for ( t = tl2; t; t = NEXT(t) ) { + /* s = [i,[*,j,*,*],...] */ + s = BDY((LIST)BDY(t)); + j = perm[QTOS((Q)ARG0(s))]; STOQ(j,jq); ARG0(s) = (pointer)jq; + for ( s = BDY((LIST)ARG1(s)); s; s = NEXT(s) ) { + j = perm[QTOS((Q)ARG1(BDY((LIST)BDY(s))))]; STOQ(j,jq); + ARG1(BDY((LIST)BDY(s))) = (pointer)jq; + } } - } - for ( j = length(x)-1, t = 0; j >= 0; j-- ) { - STOQ(perm[j],jq); MKNODE(s,jq,t); t = s; - } - MKLIST(l1,tl1); MKLIST(l2,tl2); MKLIST(l3,t); MKLIST(l4,tl3); - MKLIST(l5,tl4); - STOQ(nd_bpe,bpe); - tr = mknode(8,*rp,(!ishomo&&homo)?ONE:0,l1,l2,l3,l4,l5,bpe); MKLIST(*rp,tr); + for ( j = length(x)-1, t = 0; j >= 0; j-- ) { + STOQ(perm[j],jq); MKNODE(s,jq,t); t = s; + } + MKLIST(l1,tl1); MKLIST(l2,tl2); MKLIST(l3,t); MKLIST(l4,tl3); + MKLIST(l5,tl4); + STOQ(nd_bpe,bpe); + tr = mknode(8,*rp,(!ishomo&&homo)?ONE:0,l1,l2,l3,l4,l5,bpe); MKLIST(*rp,tr); + } } #if 0 fprintf(asir_out,"ndv_alloc=%d\n",ndv_alloc); @@ -3074,7 +3154,7 @@ void nd_gr_postproc(LIST f,LIST v,int m,struct order_s BDY(r) = ndvtop(m,CO,vv,BDY(t)); } if ( r0 ) NEXT(r) = 0; - if ( nalg ) + if ( !m && nd_nalg ) r0 = postprocess_algcoef(av,alist,r0); MKLIST(*rp,r0); } @@ -3403,7 +3483,7 @@ void nd_gr_trace(LIST f,LIST v,int trace,int homo,int if ( nd_module ) BDY(r) = ndvtopl(0,CO,vv,BDY(r),mrank); else BDY(r) = (pointer)ndvtop(0,CO,vv,BDY(r)); } - if ( nalg ) + if ( nd_nalg ) cand = postprocess_algcoef(av,alist,cand); MKLIST(*rp,cand); if ( nd_gentrace ) { @@ -4889,6 +4969,27 @@ DP ndvtodp(int mod,NDV p) return d; } +DP ndtodp(int mod,ND p) +{ + MP m,m0; + DP d; + NM t; + int i,len; + + if ( !p ) return 0; + m0 = 0; + len = p->len; + for ( t = BDY(p); t; t = NEXT(t) ) { + NEXTMP(m0,m); + m->dl = ndltodl(nd_nvar,DL(t)); + m->c = ndctop(mod,t->c); + } + NEXT(m) = 0; + MKDP(nd_nvar,m0,d); + SG(d) = SG(p); + return d; +} + void ndv_print(NDV p) { NMV m; @@ -5262,7 +5363,7 @@ Q *nm_ind_pair_to_vect(int mod,UINT *s0,int n,NM_ind_p return r; } -IndArray nm_ind_pair_to_vect_compress(int mod,UINT *s0,int n,NM_ind_pair pair) +IndArray nm_ind_pair_to_vect_compress(int mod,UINT *s0,int n,int *s0hash,NM_ind_pair pair) { NM m; NMV mr; @@ -5271,7 +5372,7 @@ IndArray nm_ind_pair_to_vect_compress(int mod,UINT *s0 unsigned char *ivc; unsigned short *ivs; UINT *v,*ivi,*s0v; - int i,j,len,prev,diff,cdiff; + int i,j,len,prev,diff,cdiff,h; IndArray r; struct oEGT eg0,eg1; @@ -5284,7 +5385,8 @@ struct oEGT eg0,eg1; get_eg(&eg0); for ( i = j = 0, s = s0, mr = BDY(p); j < len; j++, NMV_ADV(mr) ) { ndl_add(d,DL(mr),t); - for ( ; !ndl_equal(t,s); s += nd_wpd, i++ ); + h = ndl_hash_value(t); + for ( ; h != s0hash[i] || !ndl_equal(t,s); s += nd_wpd, i++ ); v[j] = i; } get_eg(&eg1); add_eg(&eg_search,&eg0,&eg1); @@ -5668,17 +5770,17 @@ int nd_symbolic_preproc(PGeoBucket bucket,int trace,UI return col; } - NODE nd_f4(int m,int **indp) { int i,nh,stat,index; - NODE r,g; - ND_pairs d,l,t; + NODE r,g,tn0,tn,node; + ND_pairs d,l,t,ll0,ll; + LIST l0,l1; ND spol,red; NDV nf,redv; NM s0,s; - NODE rp0,srp0,nflist; - int nsp,nred,col,rank,len,k,j,a; + NODE rp0,srp0,nflist,nzlist; + int nsp,nred,col,rank,len,k,j,a,i1s,i2s; UINT c; UINT **spmat; UINT *s0vect,*svect,*p,*v; @@ -5689,7 +5791,7 @@ NODE nd_f4(int m,int **indp) int sugar; PGeoBucket bucket; struct oEGT eg0,eg1,eg_f4; - + Q i1,i2,sugarq; #if 0 ndv_alloc = 0; #endif @@ -5698,10 +5800,32 @@ NODE nd_f4(int m,int **indp) d = update_pairs(d,g,i,0); g = update_base(g,i); } + nzlist = 0; while ( d ) { get_eg(&eg0); l = nd_minsugarp(d,&d); sugar = SG(l); + if ( nd_nzlist ) { + for ( tn = nd_nzlist; tn; tn = NEXT(tn) ) { + node = BDY((LIST)BDY(tn)); + if ( QTOS((Q)ARG0(node)) == sugar ) break; + } + if ( !tn ) error("nd_f4 : inconsistent non-zero list"); + for ( t = l, ll0 = 0; t; t = NEXT(t) ) { + for ( tn = BDY((LIST)ARG1(node)); tn; tn = NEXT(tn) ) { + i1s = QTOS((Q)ARG0(BDY((LIST)BDY(tn)))); + i2s = QTOS((Q)ARG1(BDY((LIST)BDY(tn)))); + if ( t->i1 == i1s && t->i2 == i2s ) break; + } + if ( tn ) { + if ( !ll0 ) ll0 = t; + else NEXT(ll) = t; + ll = t; + } + } + if ( ll0 ) NEXT(ll) = 0; + l = ll0; + } bucket = create_pbucket(); stat = nd_sp_f4(m,0,l,bucket); if ( !stat ) { @@ -5722,10 +5846,7 @@ NODE nd_f4(int m,int **indp) if ( DP_Print ) fprintf(asir_out,"sugar=%d,symb=%fsec,", sugar,eg_f4.exectime+eg_f4.gctime); - if ( 1 ) - nflist = nd_f4_red(m,l,0,s0vect,col,rp0,0); - else - nflist = nd_f4_red_dist(m,l,s0vect,col,rp0,0); + nflist = nd_f4_red(m,l,0,s0vect,col,rp0,nd_gentrace?&ll:0); /* adding new bases */ for ( r = nflist; r; r = NEXT(r) ) { nf = (NDV)BDY(r); @@ -5742,7 +5863,22 @@ NODE nd_f4(int m,int **indp) d = update_pairs(d,g,nh,0); g = update_base(g,nh); } + if ( nd_gentrace ) { + for ( t = ll, tn0 = 0; t; t = NEXT(t) ) { + NEXTNODE(tn0,tn); + STOQ(t->i1,i1); STOQ(t->i2,i2); + node = mknode(2,i1,i2); MKLIST(l0,node); + BDY(tn) = l0; + } + if ( tn0 ) NEXT(tn) = 0; MKLIST(l0,tn0); + STOQ(sugar,sugarq); node = mknode(2,sugarq,l0); MKLIST(l1,node); + MKNODE(node,l1,nzlist); nzlist = node; + } } + if ( nd_gentrace ) { + MKLIST(l0,reverse_node(nzlist)); + MKNODE(nd_alltracelist,l0,0); + } #if 0 fprintf(asir_out,"ndv_alloc=%d\n",ndv_alloc); #endif @@ -5961,6 +6097,9 @@ NODE nd_f4_red(int m,ND_pairs sp0,int trace,UINT *s0ve NODE r0,rp; ND_pairs sp; NM_ind_pair *rvect; + UINT *s; + int *s0hash; + init_eg(&eg_search); for ( sp = sp0, nsp = 0; sp; sp = NEXT(sp), nsp++ ); nred = length(rp0); @@ -5970,9 +6109,12 @@ init_eg(&eg_search); /* construction of index arrays */ rvect = (NM_ind_pair *)ALLOCA(nred*sizeof(NM_ind_pair)); + s0hash = (int *)ALLOCA(col*sizeof(int)); + for ( i = 0, s = s0vect; i < col; i++, s += nd_wpd ) + s0hash[i] = ndl_hash_value(s); for ( rp = rp0, i = 0; rp; i++, rp = NEXT(rp) ) { rvect[i] = (NM_ind_pair)BDY(rp); - imat[i] = nm_ind_pair_to_vect_compress(m,s0vect,col,rvect[i]); + imat[i] = nm_ind_pair_to_vect_compress(m,s0vect,col,s0hash,rvect[i]); rhead[imat[i]->head] = 1; } if ( m ) @@ -6042,10 +6184,17 @@ NODE nd_f4_red_main(int m,ND_pairs sp0,int nsp,UINT *s rank = nd_gauss_elim_mod(spmat,spsugar,spactive,sprow,spcol,m,colstat); r0 = 0; for ( i = 0; i < rank; i++ ) { +#if 0 NEXTNODE(r0,r); BDY(r) = (pointer)vect_to_ndv(spmat[i],spcol,col,rhead,s0vect); SG((NDV)BDY(r)) = spsugar[i]; GCFREE(spmat[i]); +#else + NEXTNODE(r0,r); BDY(r) = + (pointer)vect_to_ndv(spmat[rank-i-1],spcol,col,rhead,s0vect); + SG((NDV)BDY(r)) = spsugar[rank-i-1]; + GCFREE(spmat[rank-i-1]); +#endif } if ( r0 ) NEXT(r) = 0; @@ -6297,6 +6446,7 @@ int ox_exec_f4_red(Q proc) return s; } +#if 0 NODE nd_f4_red_dist(int m,ND_pairs sp0,UINT *s0vect,int col,NODE rp0,ND_pairs *nz) { int nsp,nred; @@ -6459,6 +6609,7 @@ void nd_exec_f4_red_dist() } fflush(nd_write); } +#endif int nd_gauss_elim_q(Q **mat0,int *sugar,int row,int col,int *colstat) { @@ -6657,6 +6808,32 @@ void ndv_save(NDV p,int index) fclose(s); } +void nd_save_mod(ND p,int index) +{ + FILE *s; + char name[BUFSIZ]; + int nv,sugar,len,c; + NM m; + + sprintf(name,"%s/%d",Demand,index); + s = fopen(name,"w"); + if ( !p ) { + len = 0; + write_int(s,&len); + fclose(s); + return; + } + nv = NV(p); + sugar = SG(p); + len = LEN(p); + write_int(s,&nv); write_int(s,&sugar); write_int(s,&len); + for ( m = BDY(p); m; m = NEXT(m) ) { + c = CM(m); write_int(s,&c); + write_intarray(s,DL(m),nd_wpd); + } + fclose(s); +} + NDV ndv_load(int index) { FILE *s; @@ -6701,6 +6878,36 @@ NDV ndv_load(int index) return d; } +ND nd_load_mod(int index) +{ + FILE *s; + char name[BUFSIZ]; + int nv,sugar,len,i,c; + ND d; + NM m0,m; + + sprintf(name,"%s/%d",Demand,index); + s = fopen(name,"r"); + /* if the file does not exist, it means p[index]=0 */ + if ( !s ) return 0; + + read_int(s,&nv); + if ( !nv ) { fclose(s); return 0; } + + read_int(s,&sugar); + read_int(s,&len); + for ( m0 = 0, i = 0; i < len; i++ ) { + NEXTNM(m0,m); + read_int(s,&c); CM(m) = c; + read_intarray(s,DL(m),nd_wpd); + } + NEXT(m) = 0; + MKND(nv,m0,len,d); + SG(d) = sugar; + fclose(s); + return d; +} + void nd_det(int mod,MAT f,P *rp) { VL fv,tv; @@ -7168,7 +7375,8 @@ void parse_nd_option(NODE opt) Obj value; nd_gentrace = 0; nd_gensyz = 0; nd_nora = 0; nd_gbblock = 0; - nd_newelim = 0; nd_intersect = 0; + nd_newelim = 0; nd_intersect = 0; nd_nzlist = 0; + nd_splist = 0; nd_check_splist = 0; for ( t = opt; t; t = NEXT(t) ) { p = BDY((LIST)BDY(t)); key = BDY((STRING)BDY(p)); @@ -7194,5 +7402,291 @@ void parse_nd_option(NODE opt) nd_newelim = value?1:0; else if ( !strcmp(key,"intersect") ) nd_intersect = value?1:0; + else if ( !strcmp(key,"trace") ) { + u = BDY((LIST)value); + nd_nzlist = BDY((LIST)ARG2(u)); + nd_bpe = QTOS((Q)ARG3(u)); + } else if ( !strcmp(key,"splist") ) + nd_splist = value?1:0; + else if ( !strcmp(key,"check_splist") ) { + nd_check_splist = BDY((LIST)value); + } } +} + +ND mdptond(DP d); +ND nd_mul_nm(int mod,NM m0,ND p); +ND *btog(NODE ti,ND **p,int nb,int mod); +ND btog_one(NODE ti,ND *p,int nb,int mod); +MAT nd_btog(LIST f,LIST v,int m,struct order_spec *ord,LIST tlist,MAT *rp); +VECT nd_btog_one(LIST f,LIST v,int m,struct order_spec *ord,LIST tlist,int pos,MAT *rp); + +/* d:monomial */ +ND mdptond(DP d) +{ + NM m; + ND r; + + if ( OID(d) == 1 ) + r = ptond(CO,CO,(P)d); + else { + NEWNM(m); + dltondl(NV(d),BDY(d)->dl,DL(m)); + CQ(m) = (Q)BDY(d)->c; + NEXT(m) = 0; + MKND(NV(d),m,1,r); + } + return r; +} + +ND nd_mul_nm(int mod,NM m0,ND p) +{ + UINT *d0; + int c0,c1,c; + NM tm,mr,mr0; + ND r; + + if ( !p ) return 0; + d0 = DL(m0); + c0 = CM(m0); + mr0 = 0; + for ( tm = BDY(p); tm; tm = NEXT(tm) ) { + NEXTNM(mr0,mr); + c = CM(tm); DMAR(c0,c,0,mod,c1); CM(mr) = c1; + ndl_add(d0,DL(tm),DL(mr)); + } + NEXT(mr) = 0; + MKND(NV(p),mr0,LEN(p),r); + return r; +} + +ND *btog(NODE ti,ND **p,int nb,int mod) +{ + PGeoBucket *r; + int i,ci; + NODE t,s; + ND m,tp; + ND *pi,*rd; + P c; + + r = (PGeoBucket *)MALLOC(nb*sizeof(PGeoBucket)); + for ( i = 0; i < nb; i++ ) + r[i] = create_pbucket(); + for ( t = ti; t; t = NEXT(t) ) { + s = BDY((LIST)BDY(t)); + if ( ARG0(s) ) { + m = mdptond((DP)ARG2(s)); + ptomp(mod,(P)HCQ(m),&c); + if ( ci = ((MQ)c)->cont ) { + HCM(m) = ci; + pi = p[QTOS((Q)ARG1(s))]; + for ( i = 0; i < nb; i++ ) { + tp = nd_mul_nm(mod,BDY(m),pi[i]); + add_pbucket(mod,r[i],tp); + } + } + ci = 1; + } else { + ptomp(mod,(P)ARG3(s),&c); ci = ((MQ)c)->cont; + ci = invm(ci,mod); + } + } + rd = (ND *)MALLOC(nb*sizeof(ND)); + for ( i = 0; i < nb; i++ ) + rd[i] = normalize_pbucket(mod,r[i]); + if ( ci != 1 ) + for ( i = 0; i < nb; i++ ) nd_mul_c(mod,rd[i],ci); + return rd; +} + +ND btog_one(NODE ti,ND *p,int nb,int mod) +{ + PGeoBucket r; + int i,ci,j; + NODE t,s; + ND m,tp; + ND pi,rd; + P c; + + r = create_pbucket(); + for ( t = ti; t; t = NEXT(t) ) { + s = BDY((LIST)BDY(t)); + if ( ARG0(s) ) { + m = mdptond((DP)ARG2(s)); + ptomp(mod,(P)HCQ(m),&c); + if ( ci = ((MQ)c)->cont ) { + HCM(m) = ci; + pi = p[j=QTOS((Q)ARG1(s))]; + if ( !pi ) { + pi = nd_load_mod(j); + tp = nd_mul_nm(mod,BDY(m),pi); + nd_free(pi); + add_pbucket(mod,r,tp); + } else { + tp = nd_mul_nm(mod,BDY(m),pi); + add_pbucket(mod,r,tp); + } + } + ci = 1; + } else { + ptomp(mod,(P)ARG3(s),&c); ci = ((MQ)c)->cont; + ci = invm(ci,mod); + } + } + rd = normalize_pbucket(mod,r); + free_pbucket(r); + if ( ci != 1 ) nd_mul_c(mod,rd,ci); + return rd; +} + +MAT nd_btog(LIST f,LIST v,int mod,struct order_spec *ord,LIST tlist,MAT *rp) +{ + int i,j,n,m,nb,pi0,pi1,nvar; + VL fv,tv,vv; + NODE permtrace,perm,trace,intred,ind,t,pi,ti; + ND **p; + ND *c; + ND u; + P inv; + MAT mat; + + parse_nd_option(current_option); + get_vars((Obj)f,&fv); pltovl(v,&vv); vlminus(fv,vv,&nd_vc); + for ( nvar = 0, tv = vv; tv; tv = NEXT(tv), nvar++ ); + switch ( ord->id ) { + case 1: + if ( ord->nv != nvar ) + error("nd_check : invalid order specification"); + break; + default: + break; + } + nd_init_ord(ord); +#if 0 + nd_bpe = QTOS((Q)ARG7(BDY(tlist))); +#else + nd_bpe = 32; +#endif + nd_setup_parameters(nvar,0); + permtrace = BDY((LIST)ARG2(BDY(tlist))); + intred = BDY((LIST)ARG3(BDY(tlist))); + ind = BDY((LIST)ARG4(BDY(tlist))); + perm = BDY((LIST)BDY(permtrace)); trace =NEXT(permtrace); + for ( i = length(perm)-1, t = trace; t; t = NEXT(t) ) { + j = QTOS((Q)BDY(BDY((LIST)BDY(t)))); + if ( j > i ) i = j; + } + n = i+1; + nb = length(BDY(f)); + p = (ND **)MALLOC(n*sizeof(ND *)); + for ( t = perm, i = 0; t; t = NEXT(t), i++ ) { + pi = BDY((LIST)BDY(t)); + pi0 = QTOS((Q)ARG0(pi)); pi1 = QTOS((Q)ARG1(pi)); + p[pi0] = c = (ND *)MALLOC(nb*sizeof(ND)); + ptomp(mod,(P)ARG2(pi),&inv); + u = ptond(CO,vv,(P)ONE); + HCM(u) = ((MQ)inv)->cont; + c[pi1] = u; + } + for ( t = trace,i=0; t; t = NEXT(t), i++ ) { + printf("%d ",i); fflush(stdout); + ti = BDY((LIST)BDY(t)); + p[j=QTOS((Q)ARG0(ti))] = btog(BDY((LIST)ARG1(ti)),p,nb,mod); + if ( j == 441 ) + printf("afo"); + } + for ( t = intred, i=0; t; t = NEXT(t), i++ ) { + printf("%d ",i); fflush(stdout); + ti = BDY((LIST)BDY(t)); + p[j=QTOS((Q)ARG0(ti))] = btog(BDY((LIST)ARG1(ti)),p,nb,mod); + if ( j == 441 ) + printf("afo"); + } + m = length(ind); + MKMAT(mat,nb,m); + for ( j = 0, t = ind; j < m; j++, t = NEXT(t) ) + for ( i = 0, c = p[QTOS((Q)BDY(t))]; i < nb; i++ ) + BDY(mat)[i][j] = ndtodp(mod,c[i]); + return mat; +} + +VECT nd_btog_one(LIST f,LIST v,int mod,struct order_spec *ord, + LIST tlist,int pos,MAT *rp) +{ + int i,j,n,m,nb,pi0,pi1,nvar; + VL fv,tv,vv; + NODE permtrace,perm,trace,intred,ind,t,pi,ti; + ND *p; + ND *c; + ND u; + P inv; + VECT vect; + + parse_nd_option(current_option); + get_vars((Obj)f,&fv); pltovl(v,&vv); vlminus(fv,vv,&nd_vc); + for ( nvar = 0, tv = vv; tv; tv = NEXT(tv), nvar++ ); + switch ( ord->id ) { + case 1: + if ( ord->nv != nvar ) + error("nd_check : invalid order specification"); + break; + default: + break; + } + nd_init_ord(ord); +#if 0 + nd_bpe = QTOS((Q)ARG7(BDY(tlist))); +#else + nd_bpe = 32; +#endif + nd_setup_parameters(nvar,0); + permtrace = BDY((LIST)ARG2(BDY(tlist))); + intred = BDY((LIST)ARG3(BDY(tlist))); + ind = BDY((LIST)ARG4(BDY(tlist))); + perm = BDY((LIST)BDY(permtrace)); trace =NEXT(permtrace); + for ( i = length(perm)-1, t = trace; t; t = NEXT(t) ) { + j = QTOS((Q)BDY(BDY((LIST)BDY(t)))); + if ( j > i ) i = j; + } + n = i+1; + nb = length(BDY(f)); + p = (ND *)MALLOC(n*sizeof(ND *)); + for ( t = perm, i = 0; t; t = NEXT(t), i++ ) { + pi = BDY((LIST)BDY(t)); + pi0 = QTOS((Q)ARG0(pi)); pi1 = QTOS((Q)ARG1(pi)); + if ( pi1 == pos ) { + ptomp(mod,(P)ARG2(pi),&inv); + u = ptond(CO,vv,(P)ONE); + HCM(u) = ((MQ)inv)->cont; + p[pi0] = u; + } + } + for ( t = trace,i=0; t; t = NEXT(t), i++ ) { + printf("%d ",i); fflush(stdout); + ti = BDY((LIST)BDY(t)); + p[j=QTOS((Q)ARG0(ti))] = btog_one(BDY((LIST)ARG1(ti)),p,nb,mod); + if ( Demand ) { + nd_save_mod(p[j],j); nd_free(p[j]); p[j] = 0; + } + } + for ( t = intred, i=0; t; t = NEXT(t), i++ ) { + printf("%d ",i); fflush(stdout); + ti = BDY((LIST)BDY(t)); + p[j=QTOS((Q)ARG0(ti))] = btog_one(BDY((LIST)ARG1(ti)),p,nb,mod); + if ( Demand ) { + nd_save_mod(p[j],j); nd_free(p[j]); p[j] = 0; + } + } + m = length(ind); + MKVECT(vect,m); + for ( j = 0, t = ind; j < m; j++, t = NEXT(t) ) { + u = p[QTOS((Q)BDY(t))]; + if ( !u ) { + u = nd_load_mod(QTOS((Q)BDY(t))); + BDY(vect)[j] = ndtodp(mod,u); + nd_free(u); + } else + BDY(vect)[j] = ndtodp(mod,u); + } + return vect; }