=================================================================== RCS file: /home/cvs/OpenXM_contrib2/asir2000/engine/nd.c,v retrieving revision 1.167 retrieving revision 1.173 diff -u -p -r1.167 -r1.173 --- OpenXM_contrib2/asir2000/engine/nd.c 2009/02/08 02:47:09 1.167 +++ OpenXM_contrib2/asir2000/engine/nd.c 2009/02/15 09:22:07 1.173 @@ -1,4 +1,4 @@ -/* $OpenXM: OpenXM_contrib2/asir2000/engine/nd.c,v 1.166 2009/02/03 08:08:01 noro Exp $ */ +/* $OpenXM: OpenXM_contrib2/asir2000/engine/nd.c,v 1.172 2009/02/15 03:07:41 noro Exp $ */ #include "nd.h" @@ -49,6 +49,7 @@ static int nd_demand; static int nd_module,nd_ispot,nd_mpos; static NODE nd_tracelist; static NODE nd_alltracelist; +static int nd_gentrace,nd_gensyz,nd_nora; NumberField get_numberfield(); UINT *nd_det_compute_bound(NDV **dm,int n,int j); @@ -67,6 +68,7 @@ P ndc_div(int mod,union oNDC a,union oNDC b); P ndctop(int mod,union oNDC c); void finalize_tracelist(int i,P cont); void conv_ilist(int demand,int trace,NODE g,int **indp); +void parse_nd_option(NODE opt); extern int Denominator,DP_Multiple; @@ -1325,7 +1327,7 @@ int nd_nf(int mod,ND d,ND g,NDV *ps,int full,NDC dn,ND p = nd_demand ? ndv_load(index) : ps[index]; /* d+g -> div*(d+g)+mul*p */ g = nd_reduce2(mod,d,g,p,mul,dn,&div); - if ( GenTrace ) { + if ( nd_gentrace ) { /* Trace=[div,index,mul,ONE] */ STOQ(index,iq); nmtodp(mod,mul,&dmul); @@ -1335,7 +1337,7 @@ int nd_nf(int mod,ND d,ND g,NDV *ps,int full,NDC dn,ND if ( !mod && g && ((double)(p_mag(HCP(g))) > hmag) ) { hg = HCU(g); nd_removecont2(d,g); - if ( dn || GenTrace ) { + if ( dn || nd_gentrace ) { /* overwrite cont : Trace=[div,index,mul,cont] */ cont = ndc_div(mod,hg,HCU(g)); if ( dn ) { @@ -1344,7 +1346,7 @@ int nd_nf(int mod,ND d,ND g,NDV *ps,int full,NDC dn,ND reductr(nd_vc,(Obj)tr,&tr1); dn->r = (R)tr1; } else divq(dn->z,(Q)cont,&dn->z); } - if ( GenTrace && !UNIQ(cont) ) ARG3(node) = (pointer)cont; + if ( nd_gentrace && !UNIQ(cont) ) ARG3(node) = (pointer)cont; } hmag = ((double)p_mag(HCP(g)))*nd_scale; } @@ -1476,36 +1478,44 @@ int nd_nf_pbucket(int mod,ND g,NDV *ps,int full,ND *rp /* input : list of NDV, cand : list of NDV */ -int ndv_check_candidate(NODE input,int obpe,int oadv,EPOS oepos,NODE cand) +int ndv_check_membership(int m,NODE input,int obpe,int oadv,EPOS oepos,NODE cand) { int n,i,stat; ND nf,d; NDV r; NODE t,s; union oNDC dn; + Q q; + LIST list; - ndv_setup(0,0,cand,0,1); + ndv_setup(m,0,cand,nd_gentrace?1:0,1); n = length(cand); + if ( nd_gentrace ) { nd_alltracelist = 0; nd_tracelist = 0; } /* membercheck : list is a subset of Id(cand) ? */ - for ( t = input; t; t = NEXT(t) ) { + for ( t = input, i = 0; t; t = NEXT(t), i++ ) { again: + nd_tracelist = 0; if ( nd_bpe > obpe ) r = ndv_dup_realloc((NDV)BDY(t),obpe,oadv,oepos); else r = (NDV)BDY(t); - d = ndvtond(0,r); - stat = nd_nf(0,0,d,nd_ps,0,0,&nf); + d = ndvtond(m,r); + stat = nd_nf(m,0,d,nd_ps,0,0,&nf); if ( !stat ) { nd_reconstruct(0,0); goto again; } else if ( nf ) return 0; + if ( nd_gentrace ) { + nd_tracelist = reverse_node(nd_tracelist); + MKLIST(list,nd_tracelist); + STOQ(i,q); s = mknode(2,q,list); MKLIST(list,s); + MKNODE(s,list,nd_alltracelist); + nd_alltracelist = s; nd_tracelist = 0; + } if ( DP_Print ) { printf("."); fflush(stdout); } } if ( DP_Print ) { printf("\n"); } - /* gbcheck : cand is a GB of Id(cand) ? */ - if ( !nd_gb(0,0,1,0) ) return 0; - /* XXX */ return 1; } @@ -1813,7 +1823,7 @@ int do_diagonalize(int sugar,int m) Q iq; for ( i = nd_psn-1; i >= 0 && SG(nd_psh[i]) == sugar; i-- ) { - if ( GenTrace ) { + if ( nd_gentrace ) { /* Trace = [1,index,1,1] */ STOQ(i,iq); node = mknode(4,ONE,iq,ONE,ONE); MKLIST(l,node); MKNODE(nd_tracelist,l,0); @@ -1829,7 +1839,7 @@ int do_diagonalize(int sugar,int m) ndv_free(nfv); hc = HCU(nf); nd_removecont(m,nf); cont = ndc_div(m,hc,HCU(nf)); - if ( GenTrace ) finalize_tracelist(i,cont); + if ( nd_gentrace ) finalize_tracelist(i,cont); nfv = ndtondv(m,nf); nd_free(nf); nd_bound[i] = ndv_compute_bound(nfv); @@ -1845,7 +1855,7 @@ int do_diagonalize(int sugar,int m) /* return value = 0 => input is not a GB */ -NODE nd_gb(int m,int ishomo,int checkonly,int **indp) +NODE nd_gb(int m,int ishomo,int checkonly,int gensyz,int **indp) { int i,nh,sugar,stat; NODE r,g,t; @@ -1861,7 +1871,7 @@ NODE nd_gb(int m,int ishomo,int checkonly,int **indp) g = 0; d = 0; for ( i = 0; i < nd_psn; i++ ) { - d = update_pairs(d,g,i); + d = update_pairs(d,g,i,gensyz); g = update_base(g,i); } sugar = 0; @@ -1888,7 +1898,7 @@ again: goto again; } #if USE_GEOBUCKET - stat = (m&&!GenTrace)?nd_nf_pbucket(m,h,nd_ps,!Top,&nf) + stat = (m&&!nd_gentrace)?nd_nf_pbucket(m,h,nd_ps,!Top,&nf) :nd_nf(m,0,h,nd_ps,!Top,0,&nf); #else stat = nd_nf(m,0,h,nd_ps,!Top,0,&nf); @@ -1898,7 +1908,7 @@ again: d = nd_reconstruct(0,d); goto again; } else if ( nf ) { - if ( checkonly ) return 0; + if ( checkonly || gensyz ) return 0; if ( DP_Print ) { printf("+"); fflush(stdout); } hc = HCU(nf); nd_removecont(m,nf); @@ -1906,7 +1916,7 @@ again: nd_monic(0,&nf); nd_removecont(m,nf); } - if ( GenTrace ) { + if ( nd_gentrace ) { cont = ndc_div(m,hc,HCU(nf)); if ( m || !UNIQ(cont) ) { t = mknode(4,0,0,0,cont); @@ -1925,10 +1935,17 @@ again: goto again; } } - d = update_pairs(d,g,nh); + d = update_pairs(d,g,nh,0); g = update_base(g,nh); FREENDP(l); } else { + if ( nd_gentrace && gensyz ) { + nd_tracelist = reverse_node(nd_tracelist); + MKLIST(list,nd_tracelist); + STOQ(-1,q); t = mknode(2,q,list); MKLIST(list,t); + MKNODE(t,list,nd_alltracelist); + nd_alltracelist = t; nd_tracelist = 0; + } if ( DP_Print ) { printf("."); fflush(stdout); } FREENDP(l); } @@ -1952,7 +1969,7 @@ int do_diagonalize_trace(int sugar,int m) P cont,cont1; for ( i = nd_psn-1; i >= 0 && SG(nd_psh[i]) == sugar; i-- ) { - if ( GenTrace ) { + if ( nd_gentrace ) { /* Trace = [1,index,1,1] */ STOQ(i,iq); node = mknode(4,ONE,iq,ONE,ONE); MKLIST(l,node); MKNODE(nd_tracelist,l,0); @@ -1979,7 +1996,7 @@ int do_diagonalize_trace(int sugar,int m) ndv_free(nfv); hc = HCU(nf); nd_removecont(0,nf); cont = ndc_div(0,hc,HCU(nf)); - if ( GenTrace ) finalize_tracelist(i,cont); + if ( nd_gentrace ) finalize_tracelist(i,cont); nfv = ndtondv(0,nf); nd_free(nf); nd_bound[i] = ndv_compute_bound(nfv); @@ -2028,7 +2045,7 @@ NODE nd_gb_trace(int m,int ishomo,int **indp) init_eg(&eg_le); g = 0; d = 0; for ( i = 0; i < nd_psn; i++ ) { - d = update_pairs(d,g,i); + d = update_pairs(d,g,i,0); g = update_base(g,i); } sugar = 0; @@ -2101,7 +2118,7 @@ again: nd_removecont(0,nfq); nfqv = ndtondv(0,nfq); nd_free(nfq); nd_removecont(m,nf); nfv = ndtondv(m,nf); nd_free(nf); } - if ( GenTrace ) { + if ( nd_gentrace ) { cont = ndc_div(0,hnfq,HCU(nfqv)); if ( !UNIQ(cont) ) { t = mknode(4,0,0,0,cont); @@ -2121,7 +2138,7 @@ again: goto again; } } - d = update_pairs(d,g,nh); + d = update_pairs(d,g,nh,0); g = update_base(g,nh); } else { if ( DP_Print ) { printf("*"); fflush(stdout); } @@ -2175,15 +2192,16 @@ NODE ndv_reduceall(int m,NODE f) union oNDC hc; P cont,cont1; + if ( nd_nora ) return f; n = length(f); ndv_setup(m,0,f,0,1); perm = (int *)MALLOC(n*sizeof(int)); - if ( GenTrace ) { + if ( nd_gentrace ) { for ( t = nd_tracelist, i = 0; i < n; i++, t = NEXT(t) ) perm[i] = QTOS((Q)ARG1(BDY((LIST)BDY(t)))); } for ( i = 0; i < n; ) { - if ( GenTrace ) { + if ( nd_gentrace ) { /* Trace = [1,index,1,1] */ STOQ(i,iq); node = mknode(4,ONE,iq,ONE,ONE); MKLIST(l,node); MKNODE(nd_tracelist,l,0); @@ -2197,7 +2215,7 @@ NODE ndv_reduceall(int m,NODE f) if ( DP_Print ) { printf("."); fflush(stdout); } ndv_free(nd_ps[i]); hc = HCU(nf); nd_removecont(m,nf); - if ( GenTrace ) { + if ( nd_gentrace ) { for ( t = nd_tracelist; t; t = NEXT(t) ) { jq = ARG1(BDY((LIST)BDY(t))); j = QTOS(jq); STOQ(perm[j],jq); ARG1(BDY((LIST)BDY(t))) = jq; @@ -2213,7 +2231,7 @@ NODE ndv_reduceall(int m,NODE f) if ( DP_Print ) { printf("\n"); } for ( a0 = 0, i = 0; i < n; i++ ) { NEXTNODE(a0,a); - if ( !GenTrace ) BDY(a) = (pointer)nd_ps[i]; + if ( !nd_gentrace ) BDY(a) = (pointer)nd_ps[i]; else { for ( j = 0; j < n; j++ ) if ( perm[j] == i ) break; BDY(a) = (pointer)nd_ps[j]; @@ -2223,16 +2241,28 @@ NODE ndv_reduceall(int m,NODE f) return a0; } -ND_pairs update_pairs( ND_pairs d, NODE /* of index */ g, int t) +ND_pairs update_pairs( ND_pairs d, NODE /* of index */ g, int t, int gensyz) { ND_pairs d1,nd,cur,head,prev,remove; if ( !g ) return d; + /* for testing */ + if ( gensyz && nd_gensyz == 2 ) { + d1 = nd_newpairs(g,t); + if ( !d ) + return d1; + else { + nd = d; + while ( NEXT(nd) ) nd = NEXT(nd); + NEXT(nd) = d1; + return d; + } + } d = crit_B(d,t); d1 = nd_newpairs(g,t); d1 = crit_M(d1); d1 = crit_F(d1); - if ( do_weyl ) + if ( gensyz || do_weyl ) head = d1; else { prev = 0; cur = head = d1; @@ -2546,7 +2576,7 @@ int ndv_newps(int m,NDV a,NDV aq) nd_ps[nd_psn] = 0; } } - if ( GenTrace ) { + if ( nd_gentrace ) { /* reverse the tracelist and append it to alltracelist */ nd_tracelist = reverse_node(nd_tracelist); MKLIST(l,nd_tracelist); STOQ(nd_psn,iq); tn = mknode(2,iq,l); MKLIST(l,tn); @@ -2617,10 +2647,10 @@ void ndv_setup(int mod,int trace,NODE f,int dont_sort, if ( mod || !dont_removecont ) ndv_removecont(mod,a); if ( !mod ) register_hcf(a); } - if ( GenTrace ) { + if ( nd_gentrace ) { STOQ(i,iq); STOQ(w[i].i,jq); node = mknode(3,iq,jq,ONE); - ARG2(node) = (pointer) - ndc_div(trace?0:mod,hc,HCU(a)); + if ( !dont_removecont ) + ARG2(node) = (pointer)ndc_div(trace?0:mod,hc,HCU(a)); MKLIST(l,node); NEXTNODE(nd_tracelist,tn); BDY(tn) = l; } NEWRHist(r); SG(r) = HTD(a); ndl_copy(HDL(a),DL(r)); @@ -2636,7 +2666,7 @@ void ndv_setup(int mod,int trace,NODE f,int dont_sort, } } } - if ( GenTrace && nd_tracelist ) NEXT(tn) = 0; + if ( nd_gentrace && nd_tracelist ) NEXT(tn) = 0; } struct order_spec *append_block(struct order_spec *spec, @@ -2758,15 +2788,18 @@ void nd_gr(LIST f,LIST v,int m,int f4,struct order_spe Obj obj; NumberField nf; struct order_spec *ord1; - NODE tr,tl1,tl2; - LIST l1,l2,l3; + NODE tr,tl1,tl2,tl3,tl4; + LIST l1,l2,l3,l4,l5; int j; Q jq; int *perm; + EPOS oepos; + int obpe,oadv,ompos; - nd_module = 0; + nd_module; if ( !m && Demand ) nd_demand = 1; else nd_demand = 0; + parse_nd_option(current_option); if ( DP_Multiple ) nd_scale = ((double)DP_Multiple)/(double)(Denominator?Denominator:1); @@ -2816,14 +2849,15 @@ void nd_gr(LIST f,LIST v,int m,int f4,struct order_spe } } nd_setup_parameters(nvar,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) ) { if ( nd_module ) { - if ( !m && !GenTrace ) pltozpl((LIST)BDY(t),&dmy,&zpl); + if ( !m && !nd_gentrace ) pltozpl((LIST)BDY(t),&dmy,&zpl); else zpl = (LIST)BDY(t); b = (pointer)pltondv(CO,vv,zpl); } else { - if ( !m && !GenTrace ) ptozp((P)BDY(t),1,&dmy,&zp); + if ( !m && !nd_gentrace ) ptozp((P)BDY(t),1,&dmy,&zp); else zp = (P)BDY(t); b = (pointer)ptondv(CO,vv,zp); } @@ -2834,14 +2868,25 @@ void nd_gr(LIST f,LIST v,int m,int f4,struct order_spe } if ( fd0 ) NEXT(fd) = 0; ndv_setup(m,0,fd0,0,0); - if ( GenTrace ) { + if ( nd_gentrace ) { MKLIST(l1,nd_tracelist); MKNODE(nd_alltracelist,l1,0); } - x = f4?nd_f4(m,&perm):nd_gb(m,ishomo,0,&perm); + x = f4?nd_f4(m,&perm):nd_gb(m,ishomo,0,0,&perm); nd_demand = 0; x = ndv_reducebase(x,perm); - if ( GenTrace ) { tl1 = nd_alltracelist; nd_alltracelist = 0; } + if ( nd_gentrace ) { tl1 = nd_alltracelist; nd_alltracelist = 0; } x = ndv_reduceall(m,x); + if ( nd_gentrace ) { + tl2 = nd_alltracelist; nd_alltracelist = 0; + ndv_check_membership(m,fd0,obpe,oadv,oepos,x); + if ( nd_gentrace ) { + tl3 = nd_alltracelist; nd_alltracelist = 0; + } else tl3 = 0; + nd_gb(m,0,1,nd_gensyz?1:0,0)!=0; + if ( nd_gentrace && nd_gensyz ) { + tl4 = nd_alltracelist; nd_alltracelist = 0; + } else tl4 = 0; + } for ( r0 = 0, t = x; t; t = NEXT(t) ) { NEXTNODE(r0,r); if ( nd_module ) BDY(r) = ndvtopl(m,CO,vv,BDY(t),mrank); @@ -2851,9 +2896,9 @@ void nd_gr(LIST f,LIST v,int m,int f4,struct order_spe if ( nalg ) r0 = postprocess_algcoef(av,alist,r0); MKLIST(*rp,r0); - if ( GenTrace ) { - tl2 = nd_alltracelist; nd_alltracelist = 0; + 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,*,*],...] */ @@ -2867,8 +2912,9 @@ void nd_gr(LIST f,LIST v,int m,int f4,struct order_spe 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); - tr = mknode(5,*rp,0,l1,l2,l3); MKLIST(*rp,tr); + MKLIST(l1,tl1); MKLIST(l2,tl2); MKLIST(l3,t); MKLIST(l4,tl3); + MKLIST(l5,tl4); + tr = mknode(7,*rp,0,l1,l2,l3,l4,l5); MKLIST(*rp,tr); } #if 0 fprintf(asir_out,"ndv_alloc=%d\n",ndv_alloc); @@ -2937,7 +2983,7 @@ void nd_gr_postproc(LIST f,LIST v,int m,struct order_s for ( x = 0, i = 0; i < nd_psn; i++ ) x = update_base(x,i); if ( do_check ) { - x = nd_gb(m,ishomo,1,&perm); + x = nd_gb(m,ishomo,1,0,&perm); if ( !x ) { *rp = 0; return; @@ -2976,13 +3022,14 @@ void nd_gr_trace(LIST f,LIST v,int trace,int homo,int NumberField nf; struct order_spec *ord1; struct oEGT eg_check,eg0,eg1; - NODE tr,tl1,tl2; - LIST l1,l2,l3; + NODE tr,tl1,tl2,tl3,tl4; + LIST l1,l2,l3,l4,l5; int *perm; - int j; + int j,ret; Q jq; nd_module = 0; + parse_nd_option(current_option); if ( DP_Multiple ) nd_scale = ((double)DP_Multiple)/(double)(Denominator?Denominator:1); @@ -3043,11 +3090,11 @@ void nd_gr_trace(LIST f,LIST v,int trace,int homo,int ishomo = 1; for ( in0 = 0, fd0 = 0, t = BDY(f); t; t = NEXT(t) ) { if ( nd_module ) { - if ( !GenTrace ) pltozpl((LIST)BDY(t),&dmy,&zpl); + if ( !nd_gentrace ) pltozpl((LIST)BDY(t),&dmy,&zpl); else zpl = (LIST)BDY(t); c = (pointer)pltondv(CO,vv,zpl); } else { - if ( !GenTrace ) ptozp((P)BDY(t),1,&dmy,&zp); + if ( !nd_gentrace ) ptozp((P)BDY(t),1,&dmy,&zp); else zp = (P)BDY(t); c = (pointer)ptondv(CO,vv,zp); } @@ -3076,7 +3123,7 @@ void nd_gr_trace(LIST f,LIST v,int trace,int homo,int if ( Demand ) nd_demand = 1; ndv_setup(m,1,fd0,0,0); - if ( GenTrace ) { + if ( nd_gentrace ) { MKLIST(l1,nd_tracelist); MKNODE(nd_alltracelist,l1,0); } cand = f4?nd_f4_trace(m,&perm):nd_gb_trace(m,ishomo || homo,&perm); @@ -3094,16 +3141,24 @@ void nd_gr_trace(LIST f,LIST v,int trace,int homo,int } nd_demand = 0; cand = ndv_reducebase(cand,perm); - if ( GenTrace ) { tl1 = nd_alltracelist; nd_alltracelist = 0; } + if ( nd_gentrace ) { tl1 = nd_alltracelist; nd_alltracelist = 0; } cand = ndv_reduceall(0,cand); cbpe = nd_bpe; - if ( GenTrace ) { tl2 = nd_alltracelist; nd_alltracelist = 0; } + if ( nd_gentrace ) { tl2 = nd_alltracelist; nd_alltracelist = 0; } + get_eg(&eg0); if ( nocheck ) break; - get_eg(&eg0); - if ( ndv_check_candidate(in0,obpe,oadv,oepos,cand) ) - /* success */ - break; + if ( ret = ndv_check_membership(0,in0,obpe,oadv,oepos,cand) ) { + if ( nd_gentrace ) { + tl3 = nd_alltracelist; nd_alltracelist = 0; + } else tl3 = 0; + /* gbcheck : cand is a GB of Id(cand) ? */ + ret = nd_gb(0,0,1,nd_gensyz?1:0,0)!=0; + if ( nd_gentrace && nd_gensyz ) { + tl4 = nd_alltracelist; nd_alltracelist = 0; + } else tl4 = 0; + } + if ( ret ) break; else if ( trace > 1 ) { /* failure */ *rp = 0; return; @@ -3133,8 +3188,9 @@ void nd_gr_trace(LIST f,LIST v,int trace,int homo,int if ( nalg ) cand = postprocess_algcoef(av,alist,cand); MKLIST(*rp,cand); - if ( GenTrace ) { + 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,*,*],...] */ @@ -3148,8 +3204,9 @@ void nd_gr_trace(LIST f,LIST v,int trace,int homo,int for ( j = length(cand)-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); - tr = mknode(5,*rp,homo?ONE:0,l1,l2,l3); MKLIST(*rp,tr); + MKLIST(l1,tl1); MKLIST(l2,tl2); MKLIST(l3,t); MKLIST(l4,tl3); + MKLIST(l5,tl4); + tr = mknode(7,*rp,(!ishomo&&homo)?ONE:0,l1,l2,l3,l4,l5); MKLIST(*rp,tr); } } @@ -3877,7 +3934,7 @@ int nd_sp(int mod,int trace,ND_pairs p,ND *rp) } t1 = ndv_mul_nm(mod,m1,p1); t2 = ndv_mul_nm(mod,m2,p2); *rp = nd_add(mod,t1,t2); - if ( GenTrace ) { + if ( nd_gentrace ) { /* nd_tracelist is initialized */ STOQ(p->i1,iq); nmtodp(mod,m1,&d); node = mknode(4,ONE,iq,d,ONE); MKLIST(hist,node); MKNODE(nd_tracelist,hist,0); @@ -5346,7 +5403,7 @@ NODE nd_f4(int m,int **indp) #endif g = 0; d = 0; for ( i = 0; i < nd_psn; i++ ) { - d = update_pairs(d,g,i); + d = update_pairs(d,g,i,0); g = update_base(g,i); } while ( d ) { @@ -5390,7 +5447,7 @@ NODE nd_f4(int m,int **indp) nf = ndtondv(m,nf1); } nh = ndv_newps(m,nf,0); - d = update_pairs(d,g,nh); + d = update_pairs(d,g,nh,0); g = update_base(g,nh); } } @@ -5424,7 +5481,7 @@ NODE nd_f4_trace(int m,int **indp) g = 0; d = 0; for ( i = 0; i < nd_psn; i++ ) { - d = update_pairs(d,g,i); + d = update_pairs(d,g,i,0); g = update_base(g,i); } while ( d ) { @@ -5490,7 +5547,7 @@ NODE nd_f4_trace(int m,int **indp) ndv_mod(m,nfv); ndv_removecont(m,nfv); nh = ndv_newps(0,nfv,nfqv); - d = update_pairs(d,g,nh); + d = update_pairs(d,g,nh,0); g = update_base(g,nh); } } @@ -6736,4 +6793,24 @@ void conv_ilist(int demand,int trace,NODE g,int **indp BDY(t) = (pointer)(demand?ndv_load(j):(trace?nd_ps_trace[j]:nd_ps[j])); } if ( indp ) *indp = ind; +} + +void parse_nd_option(NODE opt) +{ + NODE t,p; + char *key; + Obj value; + + nd_gentrace = 0; nd_gensyz = 0; nd_nora = 0; + for ( t = opt; t; t = NEXT(t) ) { + p = BDY((LIST)BDY(t)); + key = BDY((STRING)BDY(p)); + value = (Obj)BDY(NEXT(p)); + if ( !strcmp(key,"gentrace") ) + nd_gentrace = value?1:0; + else if ( !strcmp(key,"gensyz") ) + nd_gensyz = value?1:0; + else if ( !strcmp(key,"nora") ) + nd_nora = value?1:0; + } }