=================================================================== RCS file: /home/cvs/OpenXM_contrib2/asir2000/engine/nd.c,v retrieving revision 1.212 retrieving revision 1.230 diff -u -p -r1.212 -r1.230 --- OpenXM_contrib2/asir2000/engine/nd.c 2013/09/27 02:35:15 1.212 +++ OpenXM_contrib2/asir2000/engine/nd.c 2016/12/02 02:12:00 1.230 @@ -1,4 +1,4 @@ -/* $OpenXM: OpenXM_contrib2/asir2000/engine/nd.c,v 1.211 2013/09/26 08:55:11 noro Exp $ */ +/* $OpenXM: OpenXM_contrib2/asir2000/engine/nd.c,v 1.229 2016/11/17 05:33:20 noro Exp $ */ #include "nd.h" @@ -8,11 +8,14 @@ int diag_period = 6; int weight_check = 1; int (*ndl_compare_function)(UINT *a1,UINT *a2); int nd_dcomp; +int nd_rref2; NM _nm_free_list; ND _nd_free_list; ND_pairs _ndp_free_list; NODE nd_hcf; +Obj nd_top_weight; + static NODE nd_subst; static VL nd_vc; static int nd_ntrans; @@ -37,7 +40,10 @@ static UINT nd_mask[32]; static UINT nd_mask0,nd_mask1; static NDV *nd_ps; +static NDV *nd_ps_gz; static NDV *nd_ps_trace; +static NDV *nd_ps_sym; +static NDV *nd_ps_trace_sym; static RHist *nd_psh; static int nd_psn,nd_pslen; static RHist *nd_red; @@ -50,6 +56,8 @@ static int nd_found,nd_create,nd_notfirst; static int nmv_adv; static int nd_demand; static int nd_module,nd_ispot,nd_mpos,nd_pot_nelim; +static int nd_module_rank,nd_poly_weight_len; +static int *nd_poly_weight,*nd_module_weight; static NODE nd_tracelist; static NODE nd_alltracelist; static int nd_gentrace,nd_gensyz,nd_nora,nd_newelim,nd_intersect; @@ -77,9 +85,24 @@ 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); +NDV ndvtondvgz(NDV p); +NDV ndvgztondv(NDV p); +ND ndtondgz(ND p); +ND ndgztond(ND p); -extern int Denominator,DP_Multiple; +void Pdp_set_weight(NODE,VECT *); +void Pox_cmo_rpc(NODE,Obj *); +extern int Denominator,DP_Multiple,MaxDeg; + +#define BLEN (8*sizeof(unsigned long)) + +typedef struct matrix { + int row,col; + unsigned long **a; +} *matrix; + + void nd_free_private_storage() { _nm_free_list = 0; @@ -100,6 +123,20 @@ void _NM_alloc() } } +matrix alloc_matrix(int row,int col) +{ + unsigned long **a; + int i,len,blen; + matrix mat; + + mat = (matrix)MALLOC(sizeof(struct matrix)); + mat->row = row; + mat->col = col; + mat->a = a = (unsigned long **)MALLOC(row*sizeof(unsigned long *)); + return mat; +} + + void _ND_alloc() { ND p; @@ -136,6 +173,8 @@ INLINE int nd_length(ND p) } } +extern int dp_negative_weight; + INLINE int ndl_reducible(UINT *d1,UINT *d2) { UINT u1,u2; @@ -143,7 +182,7 @@ INLINE int ndl_reducible(UINT *d1,UINT *d2) if ( nd_module && (MPOS(d1) != MPOS(d2)) ) return 0; - if ( TD(d1) < TD(d2) ) return 0; + if ( !dp_negative_weight && TD(d1) < TD(d2) ) return 0; #if USE_UNROLL switch ( nd_bpe ) { case 3: @@ -490,11 +529,37 @@ int ndl_block_compare(UINT *d1,UINT *d2) int ndl_matrix_compare(UINT *d1,UINT *d2) { - int i,j,s; + int i,j,s,row; int *v; + Q **mat; + Q *w; + Q t,t1,t2; for ( j = 0; j < nd_nvar; j++ ) nd_work_vector[j] = GET_EXP(d1,j)-GET_EXP(d2,j); + if ( nd_top_weight ) { + if ( OID(nd_top_weight) == O_VECT ) { + mat = (Q **)&BDY((VECT)nd_top_weight); + row = 1; + } else { + mat = (Q **)BDY((MAT)nd_top_weight); + row = ((MAT)nd_top_weight)->row; + } + for ( i = 0; i < row; i++ ) { + w = (Q *)mat[i]; + for ( j = 0, t = 0; j < nd_nvar; j++ ) { + STOQ(nd_work_vector[j],t1); + mulq(w[j],t1,&t2); + addq(t,t2,&t1); + t = t1; + } + if ( t ) { + s = SGN(t); + if ( s > 0 ) return 1; + else if ( s < 0 ) return -1; + } + } + } for ( i = 0; i < nd_matrix_len; i++ ) { v = nd_matrix[i]; for ( j = 0, s = 0; j < nd_nvar; j++ ) @@ -502,6 +567,8 @@ int ndl_matrix_compare(UINT *d1,UINT *d2) if ( s > 0 ) return 1; else if ( s < 0 ) return -1; } + if ( !ndl_equal(d1,d2) ) + error("afo"); return 0; } @@ -584,10 +651,28 @@ int ndl_ww_lex_compare(UINT *d1,UINT *d2) return ndl_lex_compare(d1,d2); } +int ndl_module_weight_compare(UINT *d1,UINT *d2) +{ + int s,j; + + if ( nd_nvar != nd_poly_weight_len ) + error("invalid module weight : the length of polynomial weight != the number of variables"); + s = 0; + for ( j = 0; j < nd_nvar; j++ ) + s += (GET_EXP(d1,j)-GET_EXP(d2,j))*nd_poly_weight[j]; + if ( MPOS(d1) >= 1 && MPOS(d2) >= 1 ) { + s += nd_module_weight[MPOS(d1)-1]-nd_module_weight[MPOS(d2)-1]; + } + if ( s > 0 ) return 1; + else if ( s < 0 ) return -1; + else return 0; +} + int ndl_module_grlex_compare(UINT *d1,UINT *d2) { int i,c; + if ( nd_module_rank && (c = ndl_module_weight_compare(d1,d2)) ) return c; if ( nd_ispot ) { if ( nd_pot_nelim && MPOS(d1)>=nd_pot_nelim+1 && MPOS(d2) >= nd_pot_nelim+1 ) { if ( TD(d1) > TD(d2) ) return 1; @@ -614,6 +699,7 @@ int ndl_module_glex_compare(UINT *d1,UINT *d2) { int i,c; + if ( nd_module_rank && (c = ndl_module_weight_compare(d1,d2)) ) return c; if ( nd_ispot ) { if ( MPOS(d1) < MPOS(d2) ) return 1; else if ( MPOS(d1) > MPOS(d2) ) return -1; @@ -632,6 +718,7 @@ int ndl_module_lex_compare(UINT *d1,UINT *d2) { int i,c; + if ( nd_module_rank && (c = ndl_module_weight_compare(d1,d2)) ) return c; if ( nd_ispot ) { if ( MPOS(d1) < MPOS(d2) ) return 1; else if ( MPOS(d1) > MPOS(d2) ) return -1; @@ -648,6 +735,7 @@ int ndl_module_block_compare(UINT *d1,UINT *d2) { int i,c; + if ( nd_module_rank && (c = ndl_module_weight_compare(d1,d2)) ) return c; if ( nd_ispot ) { if ( MPOS(d1) < MPOS(d2) ) return 1; else if ( MPOS(d1) > MPOS(d2) ) return -1; @@ -664,6 +752,7 @@ int ndl_module_matrix_compare(UINT *d1,UINT *d2) { int i,c; + if ( nd_module_rank && (c = ndl_module_weight_compare(d1,d2)) ) return c; if ( nd_ispot ) { if ( MPOS(d1) < MPOS(d2) ) return 1; else if ( MPOS(d1) > MPOS(d2) ) return -1; @@ -680,6 +769,7 @@ int ndl_module_composite_compare(UINT *d1,UINT *d2) { int i,c; + if ( nd_module_rank && (c = ndl_module_weight_compare(d1,d2)) ) return c; if ( nd_ispot ) { if ( MPOS(d1) > MPOS(d2) ) return 1; else if ( MPOS(d1) < MPOS(d2) ) return -1; @@ -1917,6 +2007,7 @@ NODE nd_gb(int m,int ishomo,int checkonly,int gensyz,i while ( d ) { again: l = nd_minp(d,&d); + if ( MaxDeg > 0 && SG(l) > MaxDeg ) break; if ( SG(l) != sugar ) { if ( ishomo ) { diag_count = 0; @@ -1969,7 +2060,7 @@ again: } } nfv = ndtondv(m,nf); nd_free(nf); - nh = ndv_newps(m,nfv,0); + nh = ndv_newps(m,nfv,0,0); if ( !m && (ishomo && ++diag_count == diag_period) ) { diag_count = 0; stat = do_diagonalize(sugar,m); @@ -2038,6 +2129,46 @@ again: return 1; } +int check_splist_f4(int m,NODE splist) +{ + UINT *s0vect; + PGeoBucket bucket; + NODE p,rp0,t; + ND_pairs d,r,l,ll; + int col,stat; + + 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 ) { + l = nd_minsugarp(d,&d); + bucket = create_pbucket(); + stat = nd_sp_f4(m,0,l,bucket); + if ( !stat ) { + for ( ll = l; NEXT(ll); ll = NEXT(ll) ); + NEXT(ll) = d; d = l; + d = nd_reconstruct(0,d); + continue; + } + if ( bucket->m < 0 ) continue; + col = nd_symbolic_preproc(bucket,0,&s0vect,&rp0); + if ( !col ) { + for ( ll = l; NEXT(ll); ll = NEXT(ll) ); + NEXT(ll) = d; d = l; + d = nd_reconstruct(0,d); + continue; + } + if ( nd_f4_red(m,l,0,s0vect,col,rp0,0) ) return 0; + } + return 1; +} + int do_diagonalize_trace(int sugar,int m) { int i,nh,stat; @@ -2135,6 +2266,7 @@ NODE nd_gb_trace(int m,int ishomo,int **indp) while ( d ) { again: l = nd_minp(d,&d); + if ( MaxDeg > 0 && SG(l) > MaxDeg ) break; if ( SG(l) != sugar ) { #if 1 if ( ishomo ) { @@ -2209,7 +2341,7 @@ again: nd_tracelist = t; } } - nh = ndv_newps(0,nfv,nfqv); + nh = ndv_newps(0,nfv,nfqv,0); if ( ishomo && ++diag_count == diag_period ) { diag_count = 0; if ( DP_Print > 2 ) fprintf(asir_out,"|"); @@ -2232,9 +2364,11 @@ again: FREENDP(l); } if ( nd_nalg ) { - print_eg("monic",&eg_monic); - print_eg("invdalg",&eg_invdalg); - print_eg("le",&eg_le); + if ( DP_Print ) { + print_eg("monic",&eg_monic); + print_eg("invdalg",&eg_invdalg); + print_eg("le",&eg_le); + } } conv_ilist(nd_demand,1,g,indp); if ( DP_Print ) { printf("nd_gb_trace done.\n"); fflush(stdout); } @@ -2632,7 +2766,7 @@ ND_pairs nd_minsugarp( ND_pairs d, ND_pairs *prest ) return dm0; } -int ndv_newps(int m,NDV a,NDV aq) +int ndv_newps(int m,NDV a,NDV aq,int f4) { int len; RHist r; @@ -2644,15 +2778,19 @@ int ndv_newps(int m,NDV a,NDV aq) if ( nd_psn == nd_pslen ) { nd_pslen *= 2; nd_ps = (NDV *)REALLOC((char *)nd_ps,nd_pslen*sizeof(NDV)); + nd_ps_gz = (NDV *)REALLOC((char *)nd_ps_gz,nd_pslen*sizeof(NDV)); nd_ps_trace = (NDV *)REALLOC((char *)nd_ps_trace,nd_pslen*sizeof(NDV)); nd_psh = (RHist *)REALLOC((char *)nd_psh,nd_pslen*sizeof(RHist)); nd_bound = (UINT **) REALLOC((char *)nd_bound,nd_pslen*sizeof(UINT *)); + nd_ps_sym = (NDV *)REALLOC((char *)nd_ps_sym,nd_pslen*sizeof(NDV)); + nd_ps_trace_sym = (NDV *)REALLOC((char *)nd_ps_trace_sym,nd_pslen*sizeof(NDV)); } NEWRHist(r); nd_psh[nd_psn] = r; nd_ps[nd_psn] = a; if ( aq ) { nd_ps_trace[nd_psn] = aq; + if ( !nd_vc ) nd_ps_gz[nd_psn] = ndvtondvgz(aq); register_hcf(aq); nd_bound[nd_psn] = ndv_compute_bound(aq); SG(r) = SG(aq); ndl_copy(HDL(aq),DL(r)); @@ -2660,13 +2798,16 @@ int ndv_newps(int m,NDV a,NDV aq) if ( !m ) register_hcf(a); nd_bound[nd_psn] = ndv_compute_bound(a); SG(r) = SG(a); ndl_copy(HDL(a),DL(r)); + if ( !m && !nd_vc ) nd_ps_gz[nd_psn] = ndvtondvgz(a); } if ( nd_demand ) { if ( aq ) { ndv_save(nd_ps_trace[nd_psn],nd_psn); + nd_ps_trace_sym[nd_psn] = ndv_symbolic(m,nd_ps_trace[nd_psn]); nd_ps_trace[nd_psn] = 0; } else { ndv_save(nd_ps[nd_psn],nd_psn); + nd_ps_sym[nd_psn] = ndv_symbolic(m,nd_ps[nd_psn]); nd_ps[nd_psn] = 0; } } @@ -2715,7 +2856,10 @@ int ndv_setup(int mod,int trace,NODE f,int dont_sort,i } nd_pslen = 2*nd_psn; nd_ps = (NDV *)MALLOC(nd_pslen*sizeof(NDV)); + nd_ps_gz = (NDV *)MALLOC(nd_pslen*sizeof(NDV)); nd_ps_trace = (NDV *)MALLOC(nd_pslen*sizeof(NDV)); + nd_ps_sym = (NDV *)MALLOC(nd_pslen*sizeof(NDV)); + nd_ps_trace_sym = (NDV *)MALLOC(nd_pslen*sizeof(NDV)); nd_psh = (RHist *)MALLOC(nd_pslen*sizeof(RHist)); nd_bound = (UINT **)MALLOC(nd_pslen*sizeof(UINT *)); nd_hcf = 0; @@ -2732,15 +2876,17 @@ int ndv_setup(int mod,int trace,NODE f,int dont_sort,i hc = HCU(w[i].p); if ( trace ) { a = nd_ps_trace[i] = ndv_dup(0,w[i].p); + if ( !nd_vc ) nd_ps_gz[i] = ndvtondvgz(a); if ( !dont_removecont) ndv_removecont(0,a); register_hcf(a); am = nd_ps[i] = ndv_dup(mod,a); ndv_mod(mod,am); - if ( DL_COMPARE(HDL(am),HDL(a)) ) - return 0; + if ( DL_COMPARE(HDL(am),HDL(a)) ) + return 0; ndv_removecont(mod,am); } else { a = nd_ps[i] = ndv_dup(mod,w[i].p); + if ( !mod && !nd_vc ) nd_ps_gz[i] = ndvtondvgz(a); if ( mod || !dont_removecont ) ndv_removecont(mod,a); if ( !mod ) register_hcf(a); } @@ -2756,9 +2902,11 @@ int ndv_setup(int mod,int trace,NODE f,int dont_sort,i if ( nd_demand ) { if ( trace ) { ndv_save(nd_ps_trace[i],i); + nd_ps_trace_sym[i] = ndv_symbolic(mod,nd_ps_trace[i]); nd_ps_trace[i] = 0; } else { ndv_save(nd_ps[i],i); + nd_ps_sym[i] = ndv_symbolic(mod,nd_ps[i]); nd_ps[i] = 0; } } @@ -2906,6 +3054,8 @@ void nd_gr(LIST f,LIST v,int m,int homo,int retdp,int ndv_alloc = 0; #endif get_vars((Obj)f,&fv); pltovl(v,&vv); vlminus(fv,vv,&nd_vc); + if ( m && nd_vc ) + error("nd_{gr,f4} : computation over Fp(X) is unsupported. Use dp_gr_mod_main()."); for ( nvar = 0, tv = vv; tv; tv = NEXT(tv), nvar++ ); switch ( ord->id ) { case 1: @@ -2989,8 +3139,13 @@ void nd_gr(LIST f,LIST v,int m,int homo,int retdp,int return; } if ( nd_check_splist ) { - if ( check_splist(m,nd_check_splist) ) *rp = (LIST)ONE; - else *rp = 0; + if ( f4 ) { + if ( check_splist_f4(m,nd_check_splist) ) *rp = (LIST)ONE; + else *rp = 0; + } else { + 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); @@ -3415,6 +3570,7 @@ void nd_gr_trace(LIST f,LIST v,int trace,int homo,int for ( t = fd0; t; t = NEXT(t) ) ndv_homogenize((NDV)BDY(t),obpe,oadv,oepos,ompos); } + if ( MaxDeg > 0 ) nocheck = 1; while ( 1 ) { tl1 = tl2 = tl3 = tl4 = 0; if ( Demand ) @@ -4117,10 +4273,16 @@ ND_pairs nd_reconstruct(int trace,ND_pairs d) prev_ndp_free_list = _ndp_free_list; _nm_free_list = 0; _ndp_free_list = 0; - for ( i = nd_psn-1; i >= 0; i-- ) ndv_realloc(nd_ps[i],obpe,oadv,oepos); + for ( i = nd_psn-1; i >= 0; i-- ) { + ndv_realloc(nd_ps[i],obpe,oadv,oepos); + ndv_realloc(nd_ps_sym[i],obpe,oadv,oepos); + ndv_realloc(nd_ps_gz[i],obpe,oadv,oepos); + } if ( trace ) - for ( i = nd_psn-1; i >= 0; i-- ) + for ( i = nd_psn-1; i >= 0; i-- ) { ndv_realloc(nd_ps_trace[i],obpe,oadv,oepos); + ndv_realloc(nd_ps_trace_sym[i],obpe,oadv,oepos); + } s0 = 0; for ( t = d; t; t = NEXT(t) ) { NEXTND_pairs(s0,s); @@ -4641,6 +4803,48 @@ NDV ndv_dup(int mod,NDV p) return d; } +NDV ndvtondvgz(NDV p) +{ + NDV r; + int len,i; + NMV t; + + r = ndv_dup(0,p); + len = LEN(p); + for ( t = BDY(r), i = 0; i < len; i++, NMV_ADV(t) ) CZ(t) = ztogz(CQ(t)); + return r; +} + +NDV ndvgztondv(NDV p) +{ + NDV r; + int len,i; + NMV t; + + r = ndv_dup(0,p); + len = LEN(p); + for ( t = BDY(r), i = 0; i < len; i++, NMV_ADV(t) ) CQ(t) = gztoz(CZ(t)); + return r; +} + +NDV ndv_symbolic(int mod,NDV p) +{ + NDV d; + NMV t,m,m0; + int i,len; + + if ( !p ) return 0; + len = LEN(p); + m0 = m = (NMV)(mod?MALLOC_ATOMIC(len*nmv_adv):MALLOC(len*nmv_adv)); + for ( t = BDY(p), i = 0; i < len; i++, NMV_ADV(t), NMV_ADV(m) ) { + ndl_copy(DL(t),DL(m)); + CQ(m) = ONE; + } + MKNDV(NV(p),m0,len,d); + SG(d) = SG(p); + return d; +} + ND nd_dup(ND p) { ND d; @@ -4658,6 +4862,28 @@ ND nd_dup(ND p) return d; } +ND ndtondgz(ND p) +{ + ND r; + NM t; + + r = nd_dup(p); + for ( t = BDY(r); t; t = NEXT(t) ) CZ(t) = ztogz(CQ(t)); + return r; +} + + +ND ndgztond(ND p) +{ + ND r; + NM t; + + r = nd_dup(p); + for ( t = BDY(r); t; t = NEXT(t) ) CQ(t) = gztoz(CZ(t)); + return r; +} + + /* XXX if p->len == 0 then it represents 0 */ void ndv_mod(int mod,NDV p) @@ -5057,7 +5283,16 @@ NODE ndv_reducebase(NODE x,int *perm) void nd_init_ord(struct order_spec *ord) { - nd_module = (ord->id >= 256); + nd_module = (ord->id >= 256); + if ( nd_module ) { + nd_dcomp = -1; + nd_ispot = ord->ispot; + nd_pot_nelim = ord->pot_nelim; + nd_poly_weight_len = ord->nv; + nd_poly_weight = ord->top_weight; + nd_module_rank = ord->module_rank; + nd_module_weight = ord->module_top_weight; + } nd_matrix = 0; nd_matrix_len = 0; switch ( ord->id ) { @@ -5113,9 +5348,6 @@ void nd_init_ord(struct order_spec *ord) /* module order */ case 256: - nd_ispot = ord->ispot; - nd_pot_nelim = ord->pot_nelim; - nd_dcomp = -1; switch ( ord->ord.simple ) { case 0: nd_isrlex = 1; @@ -5135,17 +5367,11 @@ void nd_init_ord(struct order_spec *ord) break; case 257: /* block order */ - nd_ispot = ord->ispot; - nd_pot_nelim = ord->pot_nelim; - nd_dcomp = -1; nd_isrlex = 0; ndl_compare_function = ndl_module_block_compare; break; case 258: /* matrix order */ - nd_ispot = ord->ispot; - nd_pot_nelim = ord->pot_nelim; - nd_dcomp = -1; nd_isrlex = 0; nd_matrix_len = ord->ord.matrix.row; nd_matrix = ord->ord.matrix.matrix; @@ -5153,9 +5379,6 @@ void nd_init_ord(struct order_spec *ord) break; case 259: /* composite order */ - nd_ispot = ord->ispot; - nd_pot_nelim = ord->pot_nelim; - nd_dcomp = -1; nd_isrlex = 0; nd_worb_len = ord->ord.composite.length; nd_worb = ord->ord.composite.w_or_b; @@ -5340,6 +5563,44 @@ int nd_to_vect_q(UINT *s0,int n,ND d,Q *r) return i; } +unsigned long *nd_to_vect_2(UINT *s0,int n,int *s0hash,ND p) +{ + NM m; + unsigned long *v; + int i,j,h,size; + UINT *s,*t; + + size = sizeof(unsigned long)*(n+BLEN-1)/BLEN; + v = (unsigned long *)MALLOC_ATOMIC_IGNORE_OFF_PAGE(size); + bzero(v,size); + for ( i = j = 0, s = s0, m = BDY(p); m; j++, m = NEXT(m) ) { + t = DL(m); + h = ndl_hash_value(t); + for ( ; h != s0hash[i] || !ndl_equal(t,s); s += nd_wpd, i++ ); + v[i/BLEN] |= 1L <<(i%BLEN); + } + return v; +} + +int nd_nm_to_vect_2(UINT *s0,int n,int *s0hash,NDV p,NM m,unsigned long *v) +{ + NMV mr; + UINT *d,*t,*s; + int i,j,len,h,head; + + d = DL(m); + len = LEN(p); + t = (UINT *)ALLOCA(nd_wpd*sizeof(UINT)); + for ( i = j = 0, s = s0, mr = BDY(p); j < len; j++, NMV_ADV(mr) ) { + ndl_add(d,DL(mr),t); + h = ndl_hash_value(t); + for ( ; h != s0hash[i] || !ndl_equal(t,s); s += nd_wpd, i++ ); + if ( j == 0 ) head = i; + v[i/BLEN] |= 1L <<(i%BLEN); + } + return head; +} + Q *nm_ind_pair_to_vect(int mod,UINT *s0,int n,NM_ind_pair pair) { NM m; @@ -5378,7 +5639,7 @@ struct oEGT eg0,eg1; m = pair->mul; d = DL(m); - p = nd_ps[pair->index]; + p = nd_demand?nd_ps_sym[pair->index]:nd_ps[pair->index]; len = LEN(p); t = (UINT *)ALLOCA(nd_wpd*sizeof(UINT)); v = (unsigned int *)ALLOCA(len*sizeof(unsigned int)); @@ -5457,7 +5718,8 @@ int ndv_reduce_vect_q(Q *svect,int trace,int col,IndAr k = ivect->head; if ( svect[k] ) { maxrs = MAX(maxrs,rp0[i]->sugar); - redv = trace?nd_ps_trace[rp0[i]->index]:nd_ps[rp0[i]->index]; + redv = nd_demand?ndv_load(rp0[i]->index) + :(trace?nd_ps_trace[rp0[i]->index]:nd_ps[rp0[i]->index]); len = LEN(redv); mr = BDY(redv); igcd_cofactor(svect[k],CQ(mr),&gcd,&cs,&cr); chsgnq(cs,&mcs); @@ -5509,6 +5771,88 @@ int ndv_reduce_vect_q(Q *svect,int trace,int col,IndAr return maxrs; } +int ndv_reduce_vect_gz(GZ *gvect,int trace,int col,IndArray *imat,NM_ind_pair *rp0,int nred) +{ + int i,j,k,l,len,pos,prev,nz; + GZ cs,mcs,c1,c2,cr,gcd,t; + IndArray ivect; + unsigned char *ivc; + unsigned short *ivs; + unsigned int *ivi; + NDV redv; + NMV mr; + NODE rp; + int maxrs; + double hmag; + struct oVECT v; + + maxrs = 0; + for ( i = 0; i < col && !gvect[i]; i++ ); + if ( i == col ) return maxrs; + hmag = (double)n_bits_gz(gvect[i])*nd_scale; + for ( i = 0; i < nred; i++ ) { + ivect = imat[i]; + k = ivect->head; + if ( gvect[k] ) { + maxrs = MAX(maxrs,rp0[i]->sugar); + redv = nd_ps_gz[rp0[i]->index]; + len = LEN(redv); mr = BDY(redv); + gcdgz(gvect[k],CZ(mr),&gcd); + divsgz(gvect[k],gcd,&cs); + divsgz(CZ(mr),gcd,&cr); + chsgngz(cs,&mcs); + if ( !UNIGZ(cr) ) { + for ( j = 0; j < col; j++ ) { + mulgz(gvect[j],cr,&c1); gvect[j] = c1; + } + } + gvect[k] = 0; prev = k; + switch ( ivect->width ) { + case 1: + ivc = ivect->index.c; + for ( j = 1, NMV_ADV(mr); j < len; j++, NMV_ADV(mr) ) { + pos = prev+ivc[j]; prev = pos; + mulgz(CZ(mr),mcs,&c2); addgz(gvect[pos],c2,&t); gvect[pos] = t; + } + break; + case 2: + ivs = ivect->index.s; + for ( j = 1, NMV_ADV(mr); j < len; j++, NMV_ADV(mr) ) { + pos = prev+ivs[j]; prev = pos; + mulgz(CZ(mr),mcs,&c2); addgz(gvect[pos],c2,&t); gvect[pos] = t; + } + break; + case 4: + ivi = ivect->index.i; + for ( j = 1, NMV_ADV(mr); j < len; j++, NMV_ADV(mr) ) { + pos = prev+ivi[j]; prev = pos; + mulgz(CZ(mr),mcs,&c2); addgz(gvect[pos],c2,&t); gvect[pos] = t; + } + break; + } + for ( j = k+1; j < col && !gvect[j]; j++ ); + if ( j == col ) break; + if ( hmag && ((double)n_bits_gz(gvect[j]) > hmag) ) { + v.len = col; v.body = (pointer)gvect; gcdvgz(&v,&gcd); +#if 1 + for ( l = 0; l < col; l++ ) { divsgz(gvect[l],gcd,&t); gvect[l] = t; } +#endif + hmag = (double)n_bits_gz(gvect[j])*nd_scale; + } + } + } + for ( j = 0; j < col && !gvect[j]; j++ ); + if ( j < col ) { + v.len = col; v.body = (pointer)gvect; gcdvgz(&v,&gcd); + for ( l = 0; l < col; l++ ) { divsgz(gvect[l],gcd,&t); gvect[l] = t; } + } + if ( DP_Print ) { + fprintf(asir_out,"-"); fflush(asir_out); + } + return maxrs; +} + + int ndv_reduce_vect(int m,UINT *svect,int col,IndArray *imat,NM_ind_pair *rp0,int nred) { int i,j,k,len,pos,prev; @@ -5535,31 +5879,39 @@ int ndv_reduce_vect(int m,UINT *svect,int col,IndArray case 1: ivc = ivect->index.c; for ( j = 1, NMV_ADV(mr); j < len; j++, NMV_ADV(mr) ) { - pos = prev+ivc[j]; c1 = CM(mr); c2 = svect[pos]; - prev = pos; - DMA(c1,c,c2,up,lo); - if ( up ) { DSAB(m,up,lo,dmy,c3); svect[pos] = c3; - } else svect[pos] = lo; + pos = prev+ivc[j]; c1 = CM(mr); prev = pos; + if ( c1 ) { + c2 = svect[pos]; + DMA(c1,c,c2,up,lo); + if ( up ) { DSAB(m,up,lo,dmy,c3); svect[pos] = c3; + } else svect[pos] = lo; + } } break; case 2: ivs = ivect->index.s; for ( j = 1, NMV_ADV(mr); j < len; j++, NMV_ADV(mr) ) { - pos = prev+ivs[j]; c1 = CM(mr); c2 = svect[pos]; + pos = prev+ivs[j]; c1 = CM(mr); prev = pos; - DMA(c1,c,c2,up,lo); - if ( up ) { DSAB(m,up,lo,dmy,c3); svect[pos] = c3; - } else svect[pos] = lo; + if ( c1 ) { + c2 = svect[pos]; + DMA(c1,c,c2,up,lo); + if ( up ) { DSAB(m,up,lo,dmy,c3); svect[pos] = c3; + } else svect[pos] = lo; + } } break; case 4: ivi = ivect->index.i; for ( j = 1, NMV_ADV(mr); j < len; j++, NMV_ADV(mr) ) { - pos = prev+ivi[j]; c1 = CM(mr); c2 = svect[pos]; + pos = prev+ivi[j]; c1 = CM(mr); prev = pos; - DMA(c1,c,c2,up,lo); - if ( up ) { DSAB(m,up,lo,dmy,c3); svect[pos] = c3; - } else svect[pos] = lo; + if ( c1 ) { + c2 = svect[pos]; + DMA(c1,c,c2,up,lo); + if ( up ) { DSAB(m,up,lo,dmy,c3); svect[pos] = c3; + } else svect[pos] = lo; + } } break; } @@ -5648,6 +6000,28 @@ NDV vect_to_ndv(UINT *vect,int spcol,int col,int *rhea } } +NDV vect_to_ndv_2(unsigned long *vect,int col,UINT *s0vect) +{ + int j,k,len; + UINT *p; + NDV r; + NMV mr0,mr; + + for ( j = 0, len = 0; j < col; j++ ) if ( vect[j/BLEN] & (1L<<(j%BLEN)) ) len++; + if ( !len ) return 0; + else { + mr0 = (NMV)MALLOC_ATOMIC_IGNORE_OFF_PAGE(nmv_adv*len); + mr = mr0; + p = s0vect; + for ( j = 0; j < col; j++, p += nd_wpd ) + if ( vect[j/BLEN] & (1L<<(j%BLEN)) ) { + ndl_copy(p,DL(mr)); CM(mr) = 1; NMV_ADV(mr); + } + MKNDV(nd_nvar,mr0,len,r); + return r; + } +} + /* for preprocessed vector */ NDV vect_to_ndv_q(Q *vect,int spcol,int col,int *rhead,UINT *s0vect) @@ -5680,6 +6054,34 @@ NDV vect_to_ndv_q(Q *vect,int spcol,int col,int *rhead } } +NDV vect_to_ndv_gz(GZ *vect,int spcol,int col,int *rhead,UINT *s0vect) +{ + int j,k,len; + UINT *p; + Q c; + NDV r; + NMV mr0,mr; + + for ( j = 0, len = 0; j < spcol; j++ ) if ( vect[j] ) len++; + if ( !len ) return 0; + else { + mr0 = (NMV)MALLOC(nmv_adv*len); +#if 0 + ndv_alloc += nmv_adv*len; +#endif + mr = mr0; + p = s0vect; + for ( j = k = 0; j < col; j++, p += nd_wpd ) + if ( !rhead[j] ) { + if ( c = vect[k++] ) { + ndl_copy(p,DL(mr)); CZ(mr) = c; NMV_ADV(mr); + } + } + MKNDV(nd_nvar,mr0,len,r); + return r; + } +} + /* for plain vector */ NDV plain_vect_to_ndv_q(Q *vect,int col,UINT *s0vect) @@ -5739,7 +6141,10 @@ int nd_symbolic_preproc(PGeoBucket bucket,int trace,UI NDV *ps; s0 = 0; rp0 = 0; col = 0; - ps = trace?nd_ps_trace:nd_ps; + if ( nd_demand ) + ps = trace?nd_ps_trace_sym:nd_ps_sym; + else + ps = trace?nd_ps_trace:nd_ps; while ( 1 ) { head = remove_head_pbucket_symbolic(bucket); if ( !head ) break; @@ -5805,26 +6210,29 @@ NODE nd_f4(int m,int **indp) get_eg(&eg0); l = nd_minsugarp(d,&d); sugar = SG(l); + if ( MaxDeg > 0 && sugar > MaxDeg ) break; 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; + if ( tn ) { + nd_nzlist = NEXT(nd_nzlist); + 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; + } else l = 0; } bucket = create_pbucket(); stat = nd_sp_f4(m,0,l,bucket); @@ -5859,7 +6267,7 @@ NODE nd_f4(int m,int **indp) nd_removecont(m,nf1); nf = ndtondv(m,nf1); } - nh = ndv_newps(m,nf,0); + nh = ndv_newps(m,nf,0,1); d = update_pairs(d,g,nh,0); g = update_base(g,nh); } @@ -5882,7 +6290,7 @@ NODE nd_f4(int m,int **indp) #if 0 fprintf(asir_out,"ndv_alloc=%d\n",ndv_alloc); #endif - conv_ilist(0,0,g,indp); + conv_ilist(nd_demand,0,g,indp); return g; } @@ -5916,6 +6324,7 @@ NODE nd_f4_trace(int m,int **indp) get_eg(&eg0); l = nd_minsugarp(d,&d); sugar = SG(l); + if ( MaxDeg > 0 && sugar > MaxDeg ) break; bucket = create_pbucket(); stat = nd_sp_f4(m,0,l,bucket); if ( !stat ) { @@ -5974,7 +6383,7 @@ NODE nd_f4_trace(int m,int **indp) nfv = ndv_dup(0,nfqv); ndv_mod(m,nfv); ndv_removecont(m,nfv); - nh = ndv_newps(0,nfv,nfqv); + nh = ndv_newps(0,nfv,nfqv,1); d = update_pairs(d,g,nh,0); g = update_base(g,nh); } @@ -5982,7 +6391,7 @@ NODE nd_f4_trace(int m,int **indp) #if 0 fprintf(asir_out,"ndv_alloc=%d\n",ndv_alloc); #endif - conv_ilist(0,1,g,indp); + conv_ilist(nd_demand,1,g,indp); return g; } @@ -6077,7 +6486,7 @@ NODE nd_f4_pseudo_trace(int m,int **indp) nfv = ndv_dup(0,nfqv); ndv_mod(m,nfv); ndv_removecont(m,nfv); - nh = ndv_newps(0,nfv,nfqv); + nh = ndv_newps(0,nfv,nfqv,1); d = update_pairs(d,g,nh,0); g = update_base(g,nh); } @@ -6089,6 +6498,176 @@ NODE nd_f4_pseudo_trace(int m,int **indp) return g; } +int rref(matrix mat,int *sugar) +{ + int row,col,i,j,k,l,s,wcol,wj; + unsigned long bj; + unsigned long **a; + unsigned long *ai,*ak,*as,*t; + int *pivot; + + row = mat->row; + col = mat->col; + a = mat->a; + wcol = (col+BLEN-1)/BLEN; + pivot = (int *)MALLOC_ATOMIC(row*sizeof(int)); + i = 0; + for ( j = 0; j < col; j++ ) { + wj = j/BLEN; bj = 1L<<(j%BLEN); + for ( k = i; k < row; k++ ) + if ( a[k][wj] & bj ) break; + if ( k == row ) continue; + pivot[i] = j; + if ( k != i ) { + t = a[i]; a[i] = a[k]; a[k] = t; + s = sugar[i]; sugar[i] = sugar[k]; sugar[k] = s; + } + ai = a[i]; + for ( k = i+1; k < row; k++ ) { + ak = a[k]; + if ( ak[wj] & bj ) { + for ( l = wj; l < wcol; l++ ) + ak[l] ^= ai[l]; + sugar[k] = MAX(sugar[k],sugar[i]); + } + } + i++; + } + for ( k = i-1; k >= 0; k-- ) { + j = pivot[k]; wj = j/BLEN; bj = 1L<<(j%BLEN); + ak = a[k]; + for ( s = 0; s < k; s++ ) { + as = a[s]; + if ( as[wj] & bj ) { + for ( l = wj; l < wcol; l++ ) + as[l] ^= ak[l]; + sugar[s] = MAX(sugar[s],sugar[k]); + } + } + } + return i; +} + +void print_matrix(matrix mat) +{ + int row,col,i,j; + unsigned long *ai; + + row = mat->row; + col = mat->col; + printf("%d x %d\n",row,col); + for ( i = 0; i < row; i++ ) { + ai = mat->a[i]; + for ( j = 0; j < col; j++ ) { + if ( ai[j/BLEN] & (1L<<(j%BLEN)) ) putchar('1'); + else putchar('0'); + } + putchar('\n'); + } +} + +NDV vect_to_ndv_2(unsigned long *vect,int col,UINT *s0vect); + +void red_by_vect_2(matrix mat,int *sugar,unsigned long *v,int rhead,int rsugar) +{ + int row,col,wcol,wj,i,j; + unsigned long bj; + unsigned long *ai; + unsigned long **a; + int len; + int *pos; + + row = mat->row; + col = mat->col; + wcol = (col+BLEN-1)/BLEN; + pos = (int *)ALLOCA(wcol*sizeof(int)); + bzero(pos,wcol*sizeof(int)); + for ( i = j = 0; i < wcol; i++ ) + if ( v[i] ) pos[j++] = i;; + len = j; + wj = rhead/BLEN; + bj = 1L<a; + for ( i = 0; i < row; i++ ) { + ai = a[i]; + if ( ai[wj]&bj ) { + for ( j = 0; j < len; j++ ) + ai[pos[j]] ^= v[pos[j]]; + sugar[i] = MAX(sugar[i],rsugar); + } + } +} + +NODE nd_f4_red_2(ND_pairs sp0,UINT *s0vect,int col,NODE rp0,ND_pairs *nz) +{ + int nsp,nred,i,i0,k,rank,row; + NODE r0,rp; + ND_pairs sp; + ND spol; + NM_ind_pair rt; + int *s0hash; + UINT *s; + int *pivot,*sugar,*head; + matrix mat; + NM m; + NODE r; + struct oEGT eg0,eg1,eg2,eg_elim1,eg_elim2; + int rhead,rsugar,size; + unsigned long *v; + + get_eg(&eg0); +init_eg(&eg_search); + for ( sp = sp0, nsp = 0; sp; sp = NEXT(sp), nsp++ ); + nred = length(rp0); + mat = alloc_matrix(nsp,col); + s0hash = (int *)ALLOCA(col*sizeof(int)); + for ( i = 0, s = s0vect; i < col; i++, s += nd_wpd ) + s0hash[i] = ndl_hash_value(s); + + sugar = (int *)ALLOCA(nsp*sizeof(int)); + for ( i = 0, sp = sp0; sp; sp = NEXT(sp) ) { + nd_sp(2,0,sp,&spol); + if ( spol ) { + mat->a[i] = nd_to_vect_2(s0vect,col,s0hash,spol); + sugar[i] = SG(spol); + i++; + } + } + mat->row = i; + if ( DP_Print ) { + fprintf(asir_out,"%dx%d,",mat->row,mat->col); fflush(asir_out); + } + size = ((col+BLEN-1)/BLEN)*sizeof(unsigned long); + v = CALLOC((col+BLEN-1)/BLEN,sizeof(unsigned long)); + for ( rp = rp0, i = 0; rp; rp = NEXT(rp), i++ ) { + rt = (NM_ind_pair)BDY(rp); + bzero(v,size); + rhead = nd_nm_to_vect_2(s0vect,col,s0hash,nd_ps[rt->index],rt->mul,v); + rsugar = SG(nd_ps[rt->index])+TD(DL(rt->mul)); + red_by_vect_2(mat,sugar,v,rhead,rsugar); + } + + get_eg(&eg1); + init_eg(&eg_elim1); add_eg(&eg_elim1,&eg0,&eg1); + rank = rref(mat,sugar); + + for ( i = 0, r0 = 0; i < rank; i++ ) { + NEXTNODE(r0,r); + BDY(r) = (pointer)vect_to_ndv_2(mat->a[i],col,s0vect); + SG((NDV)BDY(r)) = sugar[i]; + } + if ( r0 ) NEXT(r) = 0; + get_eg(&eg2); + init_eg(&eg_elim2); add_eg(&eg_elim2,&eg1,&eg2); + if ( DP_Print ) { + fprintf(asir_out,"elim1=%fsec,elim2=%fsec\n", + eg_elim1.exectime+eg_elim1.gctime,eg_elim2.exectime+eg_elim2.gctime); + fflush(asir_out); + } + return r0; +} + + NODE nd_f4_red(int m,ND_pairs sp0,int trace,UINT *s0vect,int col,NODE rp0,ND_pairs *nz) { IndArray *imat; @@ -6100,6 +6679,9 @@ NODE nd_f4_red(int m,ND_pairs sp0,int trace,UINT *s0ve UINT *s; int *s0hash; + if ( m == 2 && nd_rref2 ) + return nd_f4_red_2(sp0,s0vect,col,rp0,nz); + init_eg(&eg_search); for ( sp = sp0, nsp = 0; sp; sp = NEXT(sp), nsp++ ); nred = length(rp0); @@ -6108,6 +6690,9 @@ init_eg(&eg_search); for ( i = 0; i < col; i++ ) rhead[i] = 0; /* construction of index arrays */ + if ( DP_Print ) { + fprintf(stderr,"%dx%d,",nsp+nred,col); + } 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 ) @@ -6120,8 +6705,8 @@ init_eg(&eg_search); if ( m ) r0 = nd_f4_red_main(m,sp0,nsp,s0vect,col,rvect,rhead,imat,nred,nz); else - r0 = nd_f4_red_q_main(sp0,nsp,trace,s0vect,col,rvect,rhead,imat,nred); -print_eg("search",&eg_search); + r0 = nd_f4_red_gz_main(sp0,nsp,trace,s0vect,col,rvect,rhead,imat,nred); + if ( DP_Print ) print_eg("search",&eg_search); return r0; } @@ -6184,17 +6769,10 @@ 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; @@ -6269,8 +6847,13 @@ NODE nd_f4_red_q_main(ND_pairs sp0,int nsp,int trace,U rank = nd_gauss_elim_q(spmat,spsugar,sprow,spcol,colstat); w = (pointer *)ALLOCA(rank*sizeof(pointer)); for ( i = 0; i < rank; i++ ) { +#if 0 w[rank-i-1] = (pointer)vect_to_ndv_q(spmat[i],spcol,col,rhead,s0vect); SG((NDV)w[rank-i-1]) = spsugar[i]; +#else + w[i] = (pointer)vect_to_ndv_q(spmat[i],spcol,col,rhead,s0vect); + SG((NDV)w[i]) = spsugar[i]; +#endif /* GCFREE(spmat[i]); */ } #if 0 @@ -6294,7 +6877,92 @@ NODE nd_f4_red_q_main(ND_pairs sp0,int nsp,int trace,U } return r0; } + +NODE nd_f4_red_gz_main(ND_pairs sp0,int nsp,int trace,UINT *s0vect,int col, + NM_ind_pair *rvect,int *rhead,IndArray *imat,int nred) +{ + int spcol,sprow,a; + int i,j,k,l,rank; + NODE r0,r; + ND_pairs sp; + ND spol; + GZ **spmat; + GZ *svect,*v; + int *colstat; + struct oEGT eg0,eg1,eg2,eg_f4,eg_f4_1,eg_f4_2; + int maxrs; + int *spsugar; + pointer *w; + + spcol = col-nred; + get_eg(&eg0); + /* elimination (1st step) */ + spmat = (GZ **)ALLOCA(nsp*sizeof(GZ *)); + svect = (GZ *)ALLOCA(col*sizeof(GZ)); + spsugar = (int *)ALLOCA(nsp*sizeof(Q)); + for ( a = sprow = 0, sp = sp0; a < nsp; a++, sp = NEXT(sp) ) { + nd_sp(0,trace,sp,&spol); + if ( !spol ) continue; + spol = ndtondgz(spol); + nd_to_vect_q(s0vect,col,spol,(Q *)svect); + maxrs = ndv_reduce_vect_gz(svect,trace,col,imat,rvect,nred); + for ( i = 0; i < col; i++ ) if ( svect[i] ) break; + if ( i < col ) { + spmat[sprow] = v = (GZ *)MALLOC(spcol*sizeof(GZ)); + for ( j = k = 0; j < col; j++ ) + if ( !rhead[j] ) v[k++] = svect[j]; + spsugar[sprow] = MAX(maxrs,SG(spol)); + sprow++; + } +/* nd_free(spol); */ + } + get_eg(&eg1); init_eg(&eg_f4_1); add_eg(&eg_f4_1,&eg0,&eg1); + if ( DP_Print ) { + fprintf(asir_out,"elim1=%fsec,",eg_f4_1.exectime+eg_f4_1.gctime); + fflush(asir_out); + } + /* free index arrays */ +/* for ( i = 0; i < nred; i++ ) GCFREE(imat[i]->index.c); */ + + /* elimination (2nd step) */ + colstat = (int *)ALLOCA(spcol*sizeof(int)); + rank = nd_gauss_elim_gz(spmat,spsugar,sprow,spcol,colstat); + w = (pointer *)ALLOCA(rank*sizeof(pointer)); + for ( i = 0; i < rank; i++ ) { +#if 0 + w[rank-i-1] = (pointer)vect_to_ndv_gz(spmat[i],spcol,col,rhead,s0vect); + w[rank-i-1] = ndvgztondv(w[rank-i-1]); + SG((NDV)w[rank-i-1]) = spsugar[i]; #else + w[i] = (pointer)vect_to_ndv_gz((Q *)spmat[i],spcol,col,rhead,s0vect); + w[i] = ndvgztondv(w[i]); + SG((NDV)w[i]) = spsugar[i]; +#endif +/* GCFREE(spmat[i]); */ + + } +#if 0 + qsort(w,rank,sizeof(NDV), + (int (*)(const void *,const void *))ndv_compare); +#endif + r0 = 0; + for ( i = 0; i < rank; i++ ) { + NEXTNODE(r0,r); BDY(r) = w[i]; + } + if ( r0 ) NEXT(r) = 0; + +/* for ( ; i < sprow; i++ ) GCFREE(spmat[i]); */ + get_eg(&eg2); init_eg(&eg_f4_2); add_eg(&eg_f4_2,&eg1,&eg2); + init_eg(&eg_f4); add_eg(&eg_f4,&eg0,&eg2); + if ( DP_Print ) { + fprintf(asir_out,"elim2=%fsec\n",eg_f4_2.exectime+eg_f4_2.gctime); + fprintf(asir_out,"nsp=%d,nred=%d,spmat=(%d,%d),rank=%d ", + nsp,nred,sprow,spcol,rank); + fprintf(asir_out,"%fsec\n",eg_f4.exectime+eg_f4.gctime); + } + return r0; +} +#else void printm(Q **mat,int row,int col) { int i,j; @@ -6632,6 +7300,27 @@ int nd_gauss_elim_q(Q **mat0,int *sugar,int row,int co return rank; } +int nd_gauss_elim_gz(GZ **mat0,int *sugar,int row,int col,int *colstat) +{ + int i,j,t,c,rank,inv; + int *ci,*ri; + GZ dn; + MAT m,nm; + + NEWMAT(m); m->row = row; m->col = col; m->body = (pointer **)mat0; + rank = gz_generic_gauss_elim(m,&nm,&dn,&ri,&ci); + for ( i = 0; i < row; i++ ) + for ( j = 0; j < col; j++ ) + mat0[i][j] = 0; + c = col-rank; + for ( i = 0; i < rank; i++ ) { + mat0[i][ri[i]] = dn; + for ( j = 0; j < c; j++ ) + mat0[i][ci[j]] = (GZ)BDY(nm)[i][j]; + } + return rank; +} + int nd_gauss_elim_mod(int **mat0,int *sugar,ND_pairs *spactive,int row,int col,int md,int *colstat) { int i,j,k,l,inv,a,rank,s; @@ -6990,7 +7679,9 @@ void nd_det(int mod,MAT f,P *rp) if ( mod ) ndv_mod(mod,d); chsgnq(ONE,&mone); for ( j = 0, sgn = 1; j < n; j++ ) { - if ( DP_Print ) fprintf(stderr,".",j); + if ( DP_Print ) { + fprintf(stderr,".",j); + } for ( i = j; i < n && !dm[i][j]; i++ ); if ( i == n ) { *rp = 0; @@ -7049,7 +7740,9 @@ void nd_det(int mod,MAT f,P *rp) } d = mjj; } - if ( DP_Print ) fprintf(stderr,"\n",k); + if ( DP_Print ) { + fprintf(stderr,"\n",k); + } if ( sgn < 0 ) if ( mod ) ndv_mul_c(mod,d,mod-1); @@ -7584,6 +8277,7 @@ MAT nd_btog(LIST f,LIST v,int mod,struct order_spec *o pi0 = QTOS((Q)ARG0(pi)); pi1 = QTOS((Q)ARG1(pi)); p[pi0] = c = (ND *)MALLOC(nb*sizeof(ND)); ptomp(mod,(P)ARG2(pi),&inv); + ((MQ)inv)->cont = invm(((MQ)inv)->cont,mod); u = ptond(CO,vv,(P)ONE); HCM(u) = ((MQ)inv)->cont; c[pi1] = u; @@ -7655,7 +8349,8 @@ VECT nd_btog_one(LIST f,LIST v,int mod,struct order_sp pi = BDY((LIST)BDY(t)); pi0 = QTOS((Q)ARG0(pi)); pi1 = QTOS((Q)ARG1(pi)); if ( pi1 == pos ) { - ptomp(mod,(P)ARG2(pi),&inv); + ptomp(mod,(P)ARG2(pi),&inv); + ((MQ)inv)->cont = invm(((MQ)inv)->cont,mod); u = ptond(CO,vv,(P)ONE); HCM(u) = ((MQ)inv)->cont; p[pi0] = u;