Main Page | Class List | File List | Class Members | File Members

ehrhart.c

Go to the documentation of this file.
00001 /***********************************************************************/
00002 /*                Ehrhart V4.20                                        */
00003 /*                copyright 1997, Doran Wilde                          */
00004 /*                copyright 1997-2000, Vincent Loechner                */
00005 /*       Permission is granted to copy, use, and distribute            */
00006 /*       for any commercial or noncommercial purpose under the terms   */
00007 /*       of the GNU General Public license, version 2, June 1991       */
00008 /*       (see file : LICENSING).                                       */
00009 /***********************************************************************/
00010 
00011 #include <stdio.h>
00012 #include <stdlib.h>
00013 #include <ctype.h>
00014 #include <string.h>
00015 #include <unistd.h>
00016 #include <assert.h>
00017 
00018 #include <polylib/polylib.h>
00019 
00020 
00021 /*! \class Ehrhart
00022 
00023 The following are mainly for debug purposes. You shouldn't need to
00024 change anything for daily usage...
00025 
00026 <p>
00027 
00028 you may define each macro independently 
00029 <ol>
00030 <li> #define EDEBUG minimal debug 
00031 <li> #define EDEBUG1 prints enumeration points
00032 <li> #define EDEBUG11 prints number of points
00033 <li> #define EDEBUG2 prints domains
00034 <li> #define EDEBUG21 prints more domains
00035 <li> #define EDEBUG3 prints systems of equations that are solved
00036 <li> #define EDEBUG4 prints message for degree reduction
00037 <li> #define EDEBUG5 prints result before simplification 
00038 <li> #define EDEBUG6 prints domains in Preprocess 
00039 <li> #define EDEBUG61 prints even more in Preprocess
00040 <li> #define EDEBUG62 prints domains in Preprocess2
00041 </ol>
00042 */
00043 
00044 
00045 /** 
00046     
00047 define this to print all constraints on the validity domains if not
00048 defined, only new constraints (not in validity domain given by the
00049 user) are printed
00050 
00051 */
00052 #define EPRINT_ALL_VALIDITY_CONSTRAINTS
00053 
00054 /* #define EDEBUG       */              /* minimal debug */
00055 /* #define EDEBUG1      */              /* prints enumeration points */
00056 /* #define EDEBUG11     */              /* prints number of points */
00057 /* #define EDEBUG2      */              /* prints domains */
00058 /* #define EDEBUG21     */              /* prints more domains */
00059 /* #define EDEBUG3      */              /* prints systems of equations that are solved */
00060 /* #define EDEBUG4      */              /* prints message for degree reduction */
00061 /* #define EDEBUG5      */              /* prints result before simplification */
00062 /* #define EDEBUG6      */              /* prints domains in Preprocess */
00063 /* #define EDEBUG61     */              /* prints even more in Preprocess */
00064 /* #define EDEBUG62     */              /* prints domains in Preprocess2 */
00065 
00066 
00067 /**
00068  Reduce the degree of resulting polynomials
00069 */
00070 #define REDUCE_DEGREE
00071 
00072 /** 
00073 define this to print one warning message per domain overflow these
00074 overflows should no longer happen since version 4.20
00075 */
00076 #define ALL_OVERFLOW_WARNINGS
00077 
00078 /******************* -----------END USER #DEFS-------- *********************/
00079 
00080 int overflow_warning_flag = 1;
00081 
00082 /*-------------------------------------------------------------------*/
00083 /* EHRHART POLYNOMIAL SYMBOLIC ALGEBRA SYSTEM                        */
00084 /*-------------------------------------------------------------------*/ 
00085 /** 
00086 
00087 EHRHART POLYNOMIAL SYMBOLIC ALGEBRA SYSTEM. The newly allocated enode
00088 can be freed with a simple free(x)
00089 
00090 @param type : enode type
00091 @param size : degree+1 for polynomial, period for periodic
00092 @param pos  : 1..nb_param, position of parameter            
00093 @return a newly allocated enode 
00094 
00095 */
00096 enode *new_enode(enode_type type,int size,int pos) {
00097   
00098   enode *res;
00099   int i;
00100   
00101   if(size == 0) {
00102     fprintf(stderr, "Allocating enode of size 0 !\n" );
00103     return NULL;
00104   }
00105   res = (enode *) malloc(sizeof(enode) + (size-1)*sizeof(evalue));
00106   res->type = type;
00107   res->size = size;
00108   res->pos = pos;
00109   for(i=0; i<size; i++) {
00110     value_init(res->arr[i].d);
00111     value_set_si(res->arr[i].d,0);
00112     res->arr[i].x.p = 0;
00113   }
00114   return res;
00115 } /* new_enode */
00116 
00117 /**
00118 releases all memory referenced by e.  (recursive)        
00119 @param e pointer to an evalue
00120 */
00121 void free_evalue_refs(evalue *e) {
00122   
00123   enode *p;
00124   int i;
00125   
00126   if (value_notzero_p(e->d)) {
00127     
00128     /* 'e' stores a constant */
00129     value_clear(e->d);
00130     value_clear(e->x.n);
00131     return; 
00132   }  
00133   value_clear(e->d);
00134   p = e->x.p;
00135   if (!p) return;       /* null pointer */
00136   for (i=0; i<p->size; i++) {
00137     free_evalue_refs(&(p->arr[i]));
00138   }
00139   free(p);
00140   return;
00141 } /* free_evalue_refs */
00142 
00143 /**
00144 
00145 @param e pointer to an evalue 
00146 @return description
00147 
00148 */
00149 enode *ecopy(enode *e) {
00150   
00151   enode *res;
00152   int i;
00153   
00154   res = new_enode(e->type,e->size,e->pos);
00155   for(i=0;i<e->size;++i) {
00156     value_assign(res->arr[i].d,e->arr[i].d);
00157     if(value_zero_p(res->arr[i].d))
00158       res->arr[i].x.p = ecopy(e->arr[i].x.p);
00159     else {
00160       value_init(res->arr[i].x.n);
00161       value_assign(res->arr[i].x.n,e->arr[i].x.n);
00162     }
00163   }
00164   return(res);
00165 } /* ecopy */
00166 
00167 /**
00168 
00169 @param DST destination file
00170 @param e pointer to evalue to be printed
00171 
00172 */
00173 void print_evalue(FILE *DST,evalue *e,char **pname) {
00174   
00175   if(value_notzero_p(e->d)) {    
00176     if(value_notone_p(e->d)) {
00177       value_print(DST,VALUE_FMT,e->x.n);
00178       fprintf(DST,"/");
00179       value_print(DST,VALUE_FMT,e->d);
00180     }  
00181     else {
00182       value_print(DST,VALUE_FMT,e->x.n);
00183     }
00184   }  
00185   else
00186     print_enode(DST,e->x.p,pname);
00187   return;
00188 } /* print_evalue */
00189 
00190 /** prints the enode  to DST  
00191 
00192 @param DST destination file 
00193 @param p pointer to enode  to be printed 
00194 @param pname array of strings, name of the parameters
00195 
00196 */
00197 void print_enode(FILE *DST,enode *p,char **pname) {
00198   
00199   int i;
00200   
00201   if (!p) {
00202     fprintf(DST, "NULL");
00203     return;
00204   }
00205   if (p->type == evector) {
00206     fprintf(DST, "{ ");
00207     for (i=0; i<p->size; i++) {
00208       print_evalue(DST, &p->arr[i], pname);
00209       if (i!=(p->size-1))
00210         fprintf(DST, ", ");
00211     }
00212     fprintf(DST, " }\n");
00213   }
00214   else if (p->type == polynomial) {
00215     fprintf(DST, "( ");
00216     for (i=p->size-1; i>=0; i--) {
00217       print_evalue(DST, &p->arr[i], pname);
00218       if (i==1) fprintf(DST, " * %s + ", pname[p->pos-1]);
00219       else if (i>1) 
00220         fprintf(DST, " * %s^%d + ", pname[p->pos-1], i);
00221     }
00222     fprintf(DST, " )\n");
00223   }
00224   else if (p->type == periodic) {
00225     fprintf(DST, "[ ");
00226     for (i=0; i<p->size; i++) {
00227       print_evalue(DST, &p->arr[i], pname);
00228       if (i!=(p->size-1)) fprintf(DST, ", ");
00229     }
00230     fprintf(DST," ]_%s", pname[p->pos-1]);
00231   }
00232   return;
00233 } /* print_enode */ 
00234 
00235 /**
00236 
00237 @param  e1 pointers to evalues
00238 @param  e2 pointers to evalues
00239 @return 1 (true) if they are equal, 0 (false) if not       
00240 
00241 */
00242 static int eequal(evalue *e1,evalue *e2) { 
00243  
00244     int i;
00245     enode *p1, *p2;
00246   
00247     if (value_ne(e1->d,e2->d))
00248         return 0;
00249   
00250     /* e1->d == e2->d */
00251     if (value_notzero_p(e1->d)) {    
00252         if (value_ne(e1->x.n,e2->x.n))
00253             return 0;
00254     
00255         /* e1->d == e2->d != 0  AND e1->n == e2->n */
00256         return 1;
00257     }
00258   
00259     /* e1->d == e2->d == 0 */
00260     p1 = e1->x.p;
00261     p2 = e2->x.p;
00262     if (p1->type != p2->type) return 0;
00263     if (p1->size != p2->size) return 0;
00264     if (p1->pos  != p2->pos) return 0;
00265     for (i=0; i<p1->size; i++)
00266         if (!eequal(&p1->arr[i], &p2->arr[i]) ) 
00267             return 0;
00268     return 1;
00269 } /* eequal */
00270 
00271 /** 
00272 
00273 @param e pointer to an evalue
00274 
00275 */
00276 void reduce_evalue (evalue *e) {
00277   
00278     enode *p;
00279     int i, j, k;
00280   
00281     if (value_notzero_p(e->d))
00282         return; /* a rational number, its already reduced */
00283     if(!(p = e->x.p))
00284         return; /* hum... an overflow probably occured */
00285   
00286     /* First reduce the components of p */
00287     for (i=0; i<p->size; i++)
00288         reduce_evalue(&p->arr[i]);
00289 
00290     if (p->type==periodic) {
00291     
00292         /* Try to reduce the period */
00293         for (i=1; i<=(p->size)/2; i++) {
00294             if ((p->size % i)==0) {
00295         
00296                 /* Can we reduce the size to i ? */
00297                 for (j=0; j<i; j++)
00298                     for (k=j+i; k<e->x.p->size; k+=i)
00299                         if (!eequal(&p->arr[j], &p->arr[k])) goto you_lose;
00300 
00301                 /* OK, lets do it */
00302                 for (j=i; j<p->size; j++) free_evalue_refs(&p->arr[j]);
00303                 p->size = i;
00304                 break;
00305 
00306             you_lose:   /* OK, lets not do it */
00307                 continue;
00308             }
00309         }
00310 
00311         /* Try to reduce its strength */
00312         if (p->size == 1) {
00313             value_clear(e->d);
00314             memcpy(e,&p->arr[0],sizeof(evalue));
00315             free(p);
00316         }
00317     }
00318     else if (p->type==polynomial) {
00319           
00320         /* Try to reduce the degree */
00321         for (i=p->size-1;i>=0;i--) {
00322             if (value_one_p(p->arr[i].d) && value_zero_p(p->arr[i].x.n))
00323         
00324                 /* Zero coefficient */
00325                 continue;
00326             else
00327                 break;
00328         }
00329         if (i==-1) p->size = 1;
00330         else if (i+1<p->size) p->size = i+1;
00331 
00332         /* Try to reduce its strength */
00333         if (p->size == 1) {
00334             value_clear(e->d);
00335             memcpy(e,&p->arr[0],sizeof(evalue));
00336             free(p);
00337         }
00338     }
00339 } /* reduce_evalue */
00340 
00341 
00342 /**
00343 
00344 multiplies two evalues and puts the result in res
00345 
00346 @param e1 pointer to an evalue
00347 @param e2 pointer to a constant evalue
00348 @param res  pointer to result evalue = e1 * e2
00349 
00350 */
00351 static void emul (evalue *e1,evalue *e2,evalue *res) {
00352   
00353     enode *p;
00354     int i;
00355     Value g;
00356 
00357     if (value_zero_p(e2->d)) {
00358         fprintf(stderr, "emul: ?expecting constant value\n");
00359         return;
00360     }  
00361     value_init(g);
00362     if (value_notzero_p(e1->d)) {
00363     
00364         value_init(res->x.n);
00365         /* Product of two rational numbers */
00366         value_multiply(res->d,e1->d,e2->d);
00367         value_multiply(res->x.n,e1->x.n,e2->x.n );
00368         Gcd(res->x.n, res->d,&g);
00369         if (value_notone_p(g)) {
00370             value_division(res->d,res->d,g);
00371             value_division(res->x.n,res->x.n,g);
00372         }
00373     }
00374     else { /* e1 is an expression */    
00375         value_set_si(res->d,0);
00376         p = e1->x.p;
00377         res->x.p = new_enode(p->type, p->size, p->pos);
00378         for (i=0; i<p->size; i++) {
00379             emul(&p->arr[i], e2, &(res->x.p->arr[i]) ); 
00380         }
00381     }
00382     value_clear(g);
00383     return;
00384 } /* emul */
00385 
00386 /**
00387 adds one evalue to evalue 'res. result = res + e1 
00388 
00389 @param e1 an evalue 
00390 @param res 
00391 
00392 */
00393 void eadd(evalue *e1,evalue *res) {
00394 
00395   int i;
00396   Value g,m1,m2;
00397  
00398   value_init(g);
00399   value_init(m1);
00400   value_init(m2);
00401   
00402   if (value_notzero_p(e1->d) && value_notzero_p(res->d)) {
00403     
00404     /* Add two rational numbers*/
00405     value_multiply(m1,e1->x.n,res->d);
00406     value_multiply(m2,res->x.n,e1->d);
00407     value_addto(res->x.n,m1,m2);
00408     value_multiply(res->d,e1->d,res->d);  
00409     Gcd(res->x.n,res->d,&g);
00410     if (value_notone_p(g)) {  
00411       value_division(res->d,res->d,g);
00412       value_division(res->x.n,res->x.n,g);
00413     }
00414     value_clear(g); value_clear(m1); value_clear(m2);
00415     return;
00416   }
00417   else if (value_notzero_p(e1->d) && value_zero_p(res->d)) {
00418     if (res->x.p->type==polynomial) {
00419       
00420       /* Add the constant to the constant term */
00421       eadd(e1, &res->x.p->arr[0]);
00422       value_clear(g); value_clear(m1); value_clear(m2);
00423       return;
00424     }
00425     else if (res->x.p->type==periodic) {
00426       
00427       /* Add the constant to all elements of periodic number */
00428       for (i=0; i<res->x.p->size; i++) { 
00429         eadd(e1, &res->x.p->arr[i]);
00430       }
00431       value_clear(g); value_clear(m1); value_clear(m2);
00432       return;
00433     }
00434     else {
00435       fprintf(stderr, "eadd: cannot add const with vector\n");
00436       value_clear(g); value_clear(m1); value_clear(m2);
00437       return;
00438     }
00439   }
00440   else if (value_zero_p(e1->d) && value_notzero_p(res->d)) {
00441     fprintf(stderr,"eadd: cannot add evalue to const\n");
00442     value_clear(g); value_clear(m1); value_clear(m2);
00443     return;
00444   }
00445   else {    /* ((e1->d==0) && (res->d==0)) */  
00446     if ((e1->x.p->type != res->x.p->type) ||
00447         (e1->x.p->pos  != res->x.p->pos )) {
00448       fprintf(stderr, "eadd: ?cannot add, incompatible types\n");
00449       value_clear(g); value_clear(m1); value_clear(m2);
00450       return;
00451     }
00452     if (e1->x.p->size == res->x.p->size) {
00453       for (i=0; i<res->x.p->size; i++) {
00454         eadd(&e1->x.p->arr[i], &res->x.p->arr[i]);
00455       }
00456       value_clear(g); value_clear(m1); value_clear(m2);
00457       return;
00458     }
00459     
00460     /* Sizes are different */
00461     if (res->x.p->type==polynomial) { 
00462       
00463       /* VIN100: if e1-size > res-size you have to copy e1 in a   */
00464       /* new enode and add res to that new node. If you do not do */
00465       /* that, you lose the the upper weight part of e1 !         */
00466       
00467       if(e1->x.p->size > res->x.p->size) {
00468         enode *tmp;
00469         tmp = ecopy(e1->x.p);
00470         for(i=0;i<res->x.p->size;++i) {
00471           eadd(&res->x.p->arr[i], &tmp->arr[i]);
00472           free_evalue_refs(&res->x.p->arr[i]);
00473         }
00474         res->x.p = tmp;
00475       }
00476       else {
00477         for (i=0; i<e1->x.p->size ; i++) {
00478           eadd(&e1->x.p->arr[i], &res->x.p->arr[i]);
00479         }
00480         value_clear(g); value_clear(m1); value_clear(m2);
00481         return;
00482       }
00483     }
00484     else if (res->x.p->type==periodic) {
00485       fprintf(stderr, "eadd: ?addition of different sized periodic nos\n");
00486       value_clear(g); value_clear(m1); value_clear(m2);
00487       return;
00488     }
00489     else { /* evector */
00490       fprintf(stderr, "eadd: ?cannot add vectors of different length\n");
00491       value_clear(g); value_clear(m1); value_clear(m2);
00492       return;
00493     }
00494   }
00495   value_clear(g); value_clear(m1); value_clear(m2);
00496   return;
00497 } /* eadd */
00498 
00499 /** 
00500 
00501 computes the inner product of two vectors. Result = result (evalue) =
00502 v1.v2 (dot product)
00503 
00504 @param  v1 an enode (vector)
00505 @param  v2 an enode (vector of constants)
00506 @param res result (evalue)
00507 
00508 */
00509 void edot(enode *v1,enode *v2,evalue *res) {
00510   
00511     int i;
00512     evalue tmp;
00513   
00514     if ((v1->type != evector) || (v2->type != evector)) {
00515         fprintf(stderr, "edot: ?expecting vectors\n");
00516         return;
00517     }
00518     if (v1->size != v2->size) {
00519         fprintf(stderr, "edot: ? vector lengths do not agree\n");
00520         return;
00521     }
00522     if (v1->size<=0) {
00523         value_set_si(res->d,1); /* set result to 0/1 */
00524         value_init(res->x.n);
00525         value_set_si(res->x.n,0);
00526         return;
00527     }
00528   
00529     /* vector v2 is expected to have only rational numbers in */
00530     /* the array.  No pointers. */  
00531     emul(&v1->arr[0],&v2->arr[0],res);
00532     for (i=1; i<v1->size; i++) {
00533         value_init(tmp.d);
00534         value_init(tmp.x.n);
00535      
00536         /* res = res + v1[i]*v2[i] */
00537         emul(&v1->arr[i],&v2->arr[i],&tmp);
00538         eadd(&tmp,res);
00539         free_evalue_refs(&tmp);
00540     }
00541     return;
00542 } /* edot */
00543 
00544 /**
00545 
00546 local recursive function used in the following ref contains the new
00547 position for each old index position
00548 
00549 @param e pointer to an evalue
00550 @param ref transformation Matrix
00551 
00552 */
00553 static void aep_evalue(evalue *e, int *ref) {
00554   
00555     enode *p;
00556     int i;
00557   
00558     if (value_notzero_p(e->d))
00559         return;         /* a rational number, its already reduced */
00560     if(!(p = e->x.p))
00561         return;         /* hum... an overflow probably occured */
00562   
00563     /* First check the components of p */
00564     for (i=0;i<p->size;i++)
00565         aep_evalue(&p->arr[i],ref);
00566   
00567     /* Then p itself */
00568     p->pos = ref[p->pos-1]+1;
00569     return;
00570 } /* aep_evalue */
00571 
00572 /** Comments */
00573 static void addeliminatedparams_evalue(evalue *e,Matrix *CT) {
00574         
00575     enode *p;
00576     int i, j;
00577     int *ref;
00578 
00579     if (value_notzero_p(e->d))
00580         return;          /* a rational number, its already reduced */
00581     if(!(p = e->x.p))
00582         return;          /* hum... an overflow probably occured */
00583   
00584     /* Compute ref */
00585     ref = (int *)malloc(sizeof(int)*(CT->NbRows-1));
00586     for(i=0;i<CT->NbRows-1;i++)
00587         for(j=0;j<CT->NbColumns;j++)
00588             if(value_notzero_p(CT->p[i][j])) {
00589                 ref[i] = j;
00590                 break;
00591             }
00592   
00593     /* Transform the references in e, using ref */
00594     aep_evalue(e,ref);
00595     free( ref );
00596     return;
00597 } /* addeliminatedparams_evalue */
00598 
00599 /** 
00600 
00601 This procedure finds an integer point contained in polyhedron D /
00602 first checks for positive values, then for negative values returns
00603 TRUE on success. Result is in min.  returns FALSE if no integer point
00604 is found
00605 
00606 <p>
00607 
00608 This is the maximum number of iterations for a given parameter to find
00609 a integer point inside the context. Kind of weird. cherche_min should
00610 
00611 <p>
00612 
00613 @param min
00614 @param D
00615 @param pos
00616 
00617 */
00618  //FIXME: needs to be rewritten !
00619 #define MAXITER 100
00620 int cherche_min(Value *min,Polyhedron *D,int pos) {
00621         
00622     Value binf, bsup;   /* upper & lower bound */
00623     Value i;
00624     int flag, maxiter;
00625   
00626     if(!D)
00627         return(1);
00628     if(pos > D->Dimension)
00629         return(1);
00630   
00631     value_init(binf); value_init(bsup);
00632     value_init(i);
00633   
00634 #ifdef EDEBUG61
00635     fprintf(stderr,"Entering Cherche min --> \n");
00636     fprintf(stderr,"LowerUpperBounds :\n");
00637     fprintf(stderr,"pos = %d\n",pos);
00638     fprintf(stderr,"current min = (");
00639     value_print(stderr,P_VALUE_FMT,min[0]);
00640     {int j;
00641     for(j=1;j<=D->Dimension ; j++) {
00642         fprintf(stderr,", ");
00643         value_print(stderr,P_VALUE_FMT,min[j]);
00644     }  
00645     }
00646     fprintf(stderr,")\n");
00647 #endif
00648         
00649     flag = lower_upper_bounds(pos,D,min,&binf,&bsup);
00650   
00651 #ifdef EDEBUG61
00652     fprintf(stderr, "flag = %d\n", flag);
00653     fprintf(stderr,"binf = ");
00654     value_print(stderr,P_VALUE_FMT,binf);
00655     fprintf(stderr,"\n");
00656     fprintf(stderr,"bsup = ");
00657     value_print(stderr,P_VALUE_FMT,bsup);
00658     fprintf(stderr,"\n");
00659 #endif
00660 
00661     if(flag&LB_INFINITY)
00662         value_set_si(binf,0);
00663   
00664     /* Loop from 0 (or binf if positive) to bsup */
00665     for(maxiter=0,(((flag&LB_INFINITY) || value_neg_p(binf)) ? 
00666             value_set_si(i,0) : value_assign(i,binf));
00667         ((flag&UB_INFINITY) || value_le(i,bsup)) && maxiter<MAXITER ;
00668         value_increment(i,i),maxiter++) {
00669     
00670         value_assign(min[pos],i);
00671         if(cherche_min(min,D->next,pos+1)) {
00672             value_clear(binf); value_clear(bsup);
00673             value_clear(i);
00674             return(1);
00675         }  
00676     }
00677   
00678     /* Descending loop from -1 (or bsup if negative) to binf */
00679     if((flag&LB_INFINITY) || value_neg_p(binf))
00680         for(maxiter=0,(((flag&UB_INFINITY) || value_pos_p(bsup))?
00681                 value_set_si(i,-1)
00682                 :value_assign(i,bsup));
00683             ((flag&LB_INFINITY) || value_ge(i,binf)) && maxiter<MAXITER  ;
00684             value_decrement(i,i),maxiter++) {
00685 
00686             value_assign(min[pos],i);
00687             if(cherche_min(min,D->next,pos+1)) {
00688                 value_clear(binf); value_clear(bsup);
00689                 value_clear(i);
00690                 return(1);
00691             }   
00692         }
00693     value_clear(binf); value_clear(bsup);
00694     value_clear(i);
00695 
00696     value_set_si(min[pos],0);
00697     return(0);       /* not found :-( */
00698 } /* cherche_min */
00699 
00700 /** 
00701 
00702 This procedure finds the smallest parallelepiped of size
00703 '<i>size[i]</i>' for every dimension i, contained in polyhedron D If
00704 this is not possible, an empty polyhedron is returned
00705 
00706 <p>
00707 
00708 <pre>
00709 written by vin100, 2000, for version 4.19
00710 modified 2002, version 5.10
00711 </pre>
00712 
00713 <p>
00714 
00715 It first finds the coordinates of the lexicographically smallest edge
00716 of the hypercube, obtained by transforming the constraints of D (by
00717 adding 'size' as many times as there are negative coeficients in each
00718 constraint), and finding the lexicographical min of this
00719 polyhedron. Then it builds the hypercube and returns it.
00720 
00721 <p>
00722 @param D
00723 @param size
00724 @param MAXRAYS
00725 
00726 */
00727 
00728 Polyhedron *Polyhedron_Preprocess(Polyhedron *D,Value *size,unsigned MAXRAYS)
00729 {
00730     Matrix *M;
00731     int i, j, d;
00732     Polyhedron *T, *S, *H, *C;
00733     Value *min,tmp;
00734   
00735     value_init(tmp);
00736     d = D->Dimension;
00737     if (MAXRAYS < 2*D->NbConstraints)
00738         MAXRAYS = 2*D->NbConstraints;
00739     M = Matrix_Alloc(MAXRAYS, D->Dimension+2);
00740     M->NbRows = D->NbConstraints;
00741   
00742     /* Original constraints */
00743     for(i=0;i<D->NbConstraints;i++)
00744         Vector_Copy(D->Constraint[i],M->p[i],(d+2));
00745 
00746 #ifdef EDEBUG6
00747     fprintf(stderr,"M for PreProcess : ");
00748     Matrix_Print(stderr,P_VALUE_FMT,M);
00749     fprintf(stderr,"\nsize == ");
00750     for( i=0 ; i<d ; i++ )
00751         value_print(stderr,P_VALUE_FMT,size[i]);
00752     fprintf(stderr,"\n");
00753 #endif
00754 
00755     /* Additionnal constraints */
00756     for(i=0;i<D->NbConstraints;i++) {
00757         if(value_zero_p(D->Constraint[i][0])) {
00758             fprintf(stderr,"Polyhedron_Preprocess: ");
00759             fprintf(stderr,"an equality was found where I did expect an inequality.\n");
00760             fprintf(stderr,"Trying to continue...\n");
00761             continue;
00762         }
00763         Vector_Copy(D->Constraint[i],M->p[M->NbRows],(d+2));
00764         for(j=1;j<=d;j++)
00765             if(value_neg_p(D->Constraint[i][j])) {
00766                 value_multiply(tmp,D->Constraint[i][j],size[j-1]);
00767                 value_addto(M->p[M->NbRows][d+1],M->p[M->NbRows][d+1],tmp);
00768             }
00769     
00770         /* If anything changed, add this new constraint */
00771         if(value_ne(M->p[M->NbRows][d+1],D->Constraint[i][d+1]))
00772             M->NbRows ++ ;
00773     }
00774   
00775 #ifdef EDEBUG6
00776     fprintf(stderr,"M used to find min : ");
00777     Matrix_Print(stderr,P_VALUE_FMT,M);
00778 #endif
00779   
00780     T = Constraints2Polyhedron(M,MAXRAYS);
00781     Matrix_Free(M);
00782     if (!T || emptyQ(T)) {
00783         if(T)
00784             Polyhedron_Free(T);
00785         value_clear(tmp);
00786         return(NULL);
00787     }
00788   
00789     /* Ok, now find the lexicographical min of T */
00790     min = (Value *) malloc(sizeof(Value) * (d+2));
00791     for(i=0;i<=d;i++) {
00792         value_init(min[i]);
00793         value_set_si(min[i],0);
00794     }  
00795     value_init(min[i]);
00796     value_set_si(min[i],1);
00797     C = Universe_Polyhedron(0);
00798     S = Polyhedron_Scan(T,C,MAXRAYS);
00799     Polyhedron_Free(C);
00800     Polyhedron_Free(T);
00801 
00802 #ifdef EDEBUG6
00803     for(i=0;i<=(d+1);i++) {
00804         value_print(stderr,P_VALUE_FMT,min[i]);
00805         fprintf(stderr," ,");
00806     }
00807     fprintf(stderr,"\n");
00808     Polyhedron_Print(stderr,P_VALUE_FMT,S);
00809     fprintf(stderr,"\n");
00810 #endif
00811     
00812     if (!cherche_min(min,S,1))
00813         {
00814             for(i=0;i<=(d+1);i++)
00815                 value_clear(min[i]);
00816             value_clear(tmp);
00817             return(NULL);
00818         }
00819     Domain_Free(S);
00820   
00821 #ifdef EDEBUG6
00822     fprintf(stderr,"min = ( ");
00823     value_print(stderr,P_VALUE_FMT,min[1]);
00824     for(i=2;i<=d;i++) {
00825         fprintf(stderr,", ");
00826         value_print(stderr,P_VALUE_FMT,min[i]);
00827     }  
00828     fprintf(stderr,")\n");
00829 #endif
00830   
00831     /* Min is the point from which we can construct the hypercube */
00832     M = Matrix_Alloc(d*2,d+2);
00833     for(i=0;i<d;i++) {
00834     
00835         /* Creates inequality  1 0..0 1 0..0 -min[i+1] */
00836         value_set_si(M->p[2*i][0],1);
00837         for(j=1;j<=d;j++)
00838             value_set_si(M->p[2*i][j],0);
00839         value_set_si(M->p[2*i][i+1],1);
00840         value_oppose(M->p[2*i][d+1],min[i+1]);
00841     
00842         /* Creates inequality  1 0..0 -1 0..0 min[i+1]+size -1 */
00843         value_set_si(M->p[2*i+1][0],1);
00844         for(j=1;j<=d;j++)
00845             value_set_si(M->p[2*i+1][j],0);
00846         value_set_si(M->p[2*i+1][i+1],-1);
00847         value_addto(M->p[2*i+1][d+1],min[i+1],size[i]);
00848         value_sub_int(M->p[2*i+1][d+1],M->p[2*i+1][d+1],1);
00849     }
00850   
00851 #ifdef EDEBUG6
00852     fprintf(stderr,"PolyhedronPreprocess: constraints H = ");
00853     Matrix_Print(stderr,P_VALUE_FMT,M);
00854 #endif
00855   
00856     H = Constraints2Polyhedron(M,MAXRAYS);
00857   
00858 #ifdef EDEBUG6
00859     Polyhedron_Print(stderr,P_VALUE_FMT,H);
00860     fprintf(stderr,"\n");
00861 #endif
00862   
00863     Matrix_Free(M);
00864     for(i=0;i<=(d+1);i++)
00865         value_clear(min[i]);
00866     free(min);
00867     value_clear(tmp);
00868     return(H);
00869 } /* Polyhedron_Preprocess */
00870 
00871 /** This procedure finds an hypercube of size 'size', containing
00872 polyhedron D increases size and lcm if necessary (and not "too big")
00873 If this is not possible, an empty polyhedron is returned
00874 
00875 <p>
00876 
00877 <pre> written by vin100, 2001, for version 4.19</pre>
00878 
00879 @param D
00880 @param size
00881 @param lcm
00882 @param MAXRAYS
00883 
00884 */
00885 Polyhedron *Polyhedron_Preprocess2(Polyhedron *D,Value *size,Value *lcm,unsigned MAXRAYS) {
00886   
00887     Matrix *c;
00888     Polyhedron *H;
00889     int i,j,r;
00890     Value n;      /* smallest/biggest value */
00891     Value s;    /* size in this dimension */
00892     Value tmp1,tmp2;
00893 
00894 #ifdef EDEBUG62
00895     int np;
00896 #endif
00897 
00898     value_init(n); value_init(s);
00899     value_init(tmp1); value_init(tmp2);
00900     c = Matrix_Alloc(D->Dimension*2,D->Dimension+2);
00901   
00902 #ifdef EDEBUG62
00903     fprintf(stderr,"\nPreProcess2 : starting\n");
00904     fprintf(stderr,"lcm = ");
00905     for( np=0 ; np<D->Dimension; np++ )
00906         value_print(stderr,VALUE_FMT,lcm[np]);
00907     fprintf(stderr,", size = ");
00908     for( np=0 ; np<D->Dimension; np++ )
00909         value_print(stderr,VALUE_FMT,size[np]);
00910     fprintf(stderr,"\n");
00911 #endif
00912   
00913     for(i=0;i<D->Dimension;i++) {
00914     
00915         /* Create constraint 1 0..0 1 0..0 -min */
00916         value_set_si(c->p[2*i][0],1);
00917         for(j=0;j<D->Dimension;j++)
00918             value_set_si(c->p[2*i][1+j],0);
00919         value_division(n,D->Ray[0][i+1],D->Ray[0][D->Dimension+1]);
00920         for(r=1;r<D->NbRays;r++) {
00921             value_division(tmp1,D->Ray[r][i+1],D->Ray[r][D->Dimension+1]);
00922             if(value_gt(n,tmp1)) {
00923         
00924                 /* New min */
00925                 value_division(n,D->Ray[r][i+1],D->Ray[r][D->Dimension+1]);
00926             }
00927         }
00928         value_set_si(c->p[2*i][i+1],1);
00929         value_oppose(c->p[2*i][D->Dimension+1],n);
00930     
00931         /* Create constraint 1 0..0 -1 0..0 max */
00932         value_set_si(c->p[2*i+1][0],1);
00933         for(j=0;j<D->Dimension;j++)
00934             value_set_si(c->p[2*i+1][1+j],0);
00935     
00936         /* n = (num+den-1)/den */
00937         value_addto(tmp1,D->Ray[0][i+1],D->Ray[0][D->Dimension+1]);
00938         value_sub_int(tmp1,tmp1,1);
00939         value_division(n,tmp1,D->Ray[0][D->Dimension+1]);
00940         for(r=1;r<D->NbRays;r++) {
00941             value_addto(tmp1,D->Ray[r][i+1],D->Ray[r][D->Dimension+1]);
00942             value_sub_int(tmp1,tmp1,1);
00943             value_division(tmp1,tmp1,D->Ray[r][D->Dimension+1]);
00944             if (value_lt(n,tmp1)) {
00945         
00946                 /* New max */
00947                 value_addto(tmp1,D->Ray[r][i+1],D->Ray[r][D->Dimension+1]);
00948                 value_sub_int(tmp1,tmp1,1);
00949                 value_division(n,tmp1,D->Ray[r][D->Dimension+1]);
00950             }
00951         }
00952         value_set_si(c->p[2*i+1][i+1],-1);
00953         value_assign(c->p[2*i+1][D->Dimension+1],n);
00954         value_addto(s,c->p[2*i+1][D->Dimension+1],c->p[2*i][D->Dimension+1]);
00955     
00956         /* Now test if the dimension of the cube is greater than the size */
00957         if(value_gt(s,size[i])) {
00958       
00959 #ifdef EDEBUG62
00960             fprintf(stderr,"size on dimension %d\n",i);
00961             fprintf(stderr,"lcm = ");
00962             for( np=0 ; np<D->Dimension; np++ )
00963                 value_print(stderr,VALUE_FMT,lcm[np]);
00964             fprintf(stderr,", size = ");
00965             for( np=0 ; np<D->Dimension; np++ )
00966                 value_print(stderr,VALUE_FMT,size[np]);
00967             fprintf(stderr,"\n");
00968             fprintf(stderr,"required size (s) = ");
00969             value_print(stderr,VALUE_FMT,s);
00970             fprintf(stderr,"\n");
00971 #endif
00972       
00973             /* If the needed size is "small enough"(<=20 or less than twice *size) */
00974             /* then increase *size, and artificially increase lcm too !*/
00975             value_set_si(tmp1,20);
00976             value_addto(tmp2,size[i],size[i]);
00977             if(value_le(s,tmp1)
00978                     || value_le(s,tmp2)) {
00979 
00980                 if( value_zero_p(lcm[i]) )
00981                     value_set_si(lcm[i],1);
00982                 /* lcm divides size... */
00983                 value_division(tmp1,size[i],lcm[i]);
00984         
00985                 /* lcm = ceil(s/h) */
00986                 value_addto(tmp2,s,tmp1);
00987                 value_add_int(tmp2,tmp2,1);
00988                 value_division(lcm[i],tmp2,tmp1);
00989         
00990                 /* new size = lcm*h */
00991                 value_multiply(size[i],lcm[i],tmp1);
00992         
00993 #ifdef EDEBUG62
00994                 fprintf(stderr,"new size = ");
00995                 for( np=0 ; np<D->Dimension; np++ )
00996                     value_print(stderr,VALUE_FMT,size[np]);
00997                 fprintf(stderr,", new lcm = ");
00998                 for( np=0 ; np<D->Dimension; np++ )
00999                     value_print(stderr,VALUE_FMT,lcm[np]);
01000                 fprintf(stderr,"\n");
01001 #endif
01002 
01003             }
01004             else {
01005         
01006 #ifdef EDEBUG62
01007                 fprintf(stderr,"Failed on dimension %d.\n",i);
01008 #endif
01009                 break;
01010             }
01011         }
01012     }
01013     if(i!=D->Dimension) {
01014         Matrix_Free(c);
01015         value_clear(n); value_clear(s);
01016         value_clear(tmp1); value_clear(tmp2);
01017         return(NULL);
01018     }
01019     for(i=0;i<D->Dimension;i++) {
01020         value_substract(c->p[2*i+1][D->Dimension+1],size[i], 
01021                 c->p[2*i][D->Dimension+1]);
01022     }
01023   
01024 #ifdef EDEBUG62
01025     fprintf(stderr,"PreProcess2 : c =");
01026     Matrix_Print(stderr,P_VALUE_FMT,c);
01027 #endif
01028 
01029     H = Constraints2Polyhedron(c,MAXRAYS);
01030     Matrix_Free(c);
01031     value_clear(n); value_clear(s);
01032     value_clear(tmp1); value_clear(tmp2);
01033     return(H);
01034 } /* Polyhedron_Preprocess2 */
01035 
01036 /**  
01037 
01038 This procedure adds additional constraints to D so that as each
01039 parameter is scanned, it will have a minimum of 'size' points If this
01040 is not possible, an empty polyhedron is returned
01041 
01042 @param D
01043 @param size
01044 @param MAXRAYS
01045 
01046 */
01047 Polyhedron *old_Polyhedron_Preprocess(Polyhedron *D,Value size,unsigned MAXRAYS) {
01048   
01049     int p, p1, ub, lb;
01050     Value a, a1, b, b1, g, aa;
01051     Value abs_a, abs_b, size_copy;
01052     int dim, con, new, needed;
01053     Value **C;
01054     Matrix *M;
01055     Polyhedron *D1;
01056   
01057     value_init(a); value_init(a1); value_init(b);
01058     value_init(b1); value_init(g); value_init(aa);
01059     value_init(abs_a); value_init(abs_b); value_init(size_copy);
01060   
01061     dim = D->Dimension;
01062     con = D->NbConstraints;
01063     M = Matrix_Alloc(MAXRAYS,dim+2);
01064     new = 0;
01065     value_assign(size_copy,size);
01066     C = D->Constraint;
01067     for (p=1; p<=dim; p++) {
01068         for (ub=0; ub<con; ub++) {
01069             value_assign(a,C[ub][p]);
01070             if (value_posz_p(a))        /* a >= 0 */
01071                 continue;               /* not an upper bound */
01072             for (lb=0;lb<con;lb++) {
01073                 value_assign(b,C[lb][p]);
01074                 if (value_negz_p(b))
01075                     continue;   /* not a lower bound */
01076         
01077                 /* Check if a new constraint is needed for this (ub,lb) pair */
01078                 /* a constraint is only needed if a previously scanned */
01079                 /* parameter (1..p-1) constrains this parameter (p) */
01080                 needed=0;
01081                 for (p1=1; p1<p; p1++) {
01082                     if (value_notzero_p(C[ub][p1]) || value_notzero_p(C[lb][p1])) {
01083                         needed=1; 
01084                         break;
01085                     }
01086                 }
01087                 if (!needed) continue;  
01088                 value_absolute(abs_a,a);
01089                 value_absolute(abs_b,b);
01090 
01091                 /* Create new constraint: b*UB-a*LB >= a*b*size */
01092                 Gcd(abs_a,abs_b,&g);
01093                 value_division(a1,a,g);
01094                 value_division(b1,b,g);
01095                 value_set_si(M->p[new][0],1);
01096                 value_oppose(abs_a,a1);           /* abs_a = -a1 */
01097                 Vector_Combine(&(C[ub][1]),&(C[lb][1]),&(M->p[new][1]),
01098                         b1,abs_a,dim+1);
01099                 value_multiply(aa,a1,b1);
01100                 value_multiply(abs_a,aa,size_copy);
01101                 value_addto(M->p[new][dim+1],M->p[new][dim+1],abs_a);
01102                 Vector_Normalize(&(M->p[new][1]),(dim+1));
01103                 new++;
01104             }
01105         }
01106     }
01107     D1 = AddConstraints(M->p_Init,new,D,MAXRAYS);
01108     Matrix_Free(M);
01109   
01110     value_clear(a); value_clear(a1); value_clear(b);
01111     value_clear(b1); value_clear(g); value_clear(aa);
01112     value_clear(abs_a); value_clear(abs_b); value_clear(size_copy);
01113     return D1;
01114 } /* old_Polyhedron_Preprocess */
01115 
01116 /** 
01117 
01118 PROCEDURES TO COMPUTE ENUMERATION. recursive procedure, recurse for
01119 each imbriquation
01120 
01121 @param pos index position of current loop index (1..hdim-1)
01122 @param P loop domain
01123 @param context context values for fixed indices 
01124 @return the number of integer points in this
01125 polyhedron 
01126 
01127 */
01128 void count_points (int pos,Polyhedron *P,Value *context, Value *res) {
01129    
01130     Value LB, UB, k, c;
01131 
01132     if (emptyQ(P)) {
01133         value_set_si(*res, 0);
01134         return;
01135     }
01136 
01137     value_init(LB); value_init(UB); value_init(k);
01138     value_set_si(LB,0);
01139     value_set_si(UB,0);
01140 
01141     if (lower_upper_bounds(pos,P,context,&LB,&UB) !=0) {
01142     
01143         /* Problem if UB or LB is INFINITY */
01144         fprintf(stderr, "count_points: ? infinite domain\n");
01145         value_clear(LB); value_clear(UB); value_clear(k);
01146         value_set_si(*res, -1);
01147         return;
01148     }
01149   
01150 #ifdef EDEBUG1
01151     if (!P->next) {
01152         int i;
01153         for (value_assign(k,LB); value_le(k,UB); value_increment(k,k)) {
01154             fprintf(stderr, "(");
01155             for (i=1; i<pos; i++) {
01156                 value_print(stderr,P_VALUE_FMT,context[i]);
01157                 fprintf(stderr,",");
01158             }
01159             value_print(stderr,P_VALUE_FMT,k);
01160             fprintf(stderr,")\n");
01161         }
01162     }
01163 #endif
01164   
01165     value_set_si(context[pos],0);
01166     if (value_lt(UB,LB)) {
01167         value_clear(LB); value_clear(UB); value_clear(k);
01168         value_set_si(*res, 0);
01169         return;  
01170     }  
01171     if (!P->next) {
01172         value_substract(k,UB,LB);
01173         value_add_int(k,k,1);
01174         value_assign(*res, k);
01175         value_clear(LB); value_clear(UB); value_clear(k);
01176         return;
01177     } 
01178 
01179     /*-----------------------------------------------------------------*/
01180     /* Optimization idea                                               */
01181     /*   If inner loops are not a function of k (the current index)    */
01182     /*   i.e. if P->Constraint[i][pos]==0 for all P following this and */
01183     /*        for all i,                                               */
01184     /*   Then CNT = (UB-LB+1)*count_points(pos+1, P->next, context)    */
01185     /*   (skip the for loop)                                           */
01186     /*-----------------------------------------------------------------*/
01187   
01188     value_init(c);
01189     value_set_si(*res, 0);
01190     for (value_assign(k,LB);value_le(k,UB);value_increment(k,k)) {
01191         /* Insert k in context */
01192         value_assign(context[pos],k);
01193         count_points(pos+1,P->next,context,&c);
01194         if(value_notmone_p(c))
01195             value_addto(*res, *res, c);
01196         else {
01197             value_set_si(*res, -1);
01198             break;
01199         }
01200     }
01201     value_clear(c);
01202   
01203 #ifdef EDEBUG11
01204     fprintf(stderr,"%d\n",CNT);
01205 #endif
01206   
01207     /* Reset context */
01208     value_set_si(context[pos],0);
01209     value_clear(LB); value_clear(UB); value_clear(k);
01210     return;
01211 } /* count_points */
01212 
01213 /*-------------------------------------------------------------------*/
01214 /* enode *P_Enum(L, LQ, context, pos, nb_param, dim, lcm,param_name) */
01215 /*     L : list of polyhedra for the loop nest                       */
01216 /*     LQ : list of polyhedra for the parameter loop nest            */
01217 /*     pos : 1..nb_param, position of the parameter                  */
01218 /*     nb_param : number of parameters total                         */
01219 /*     dim : total dimension of the polyhedron, param incl.          */
01220 /*     lcm : denominator array [0..dim-1] of the polyhedron          */
01221 /*     param_name : name of the parameters                           */
01222 /* Returns an enode tree representing the pseudo polynomial          */
01223 /*   expression for the enumeration of the polyhedron.               */
01224 /* A recursive procedure.                                            */
01225 /*-------------------------------------------------------------------*/
01226 static enode *P_Enum(Polyhedron *L,Polyhedron *LQ,Value *context,int pos,int nb_param,int dim,Value *lcm,char **param_name) {
01227 
01228   enode *res,*B,*C;
01229   int hdim,i,j,rank,flag;
01230   Value n,g,nLB,nUB,nlcm,noff,nexp,k1,nm,hdv,k,lcm_copy;
01231   Value tmp;
01232   Matrix *A;
01233 
01234 #ifdef EDEBUG
01235   fprintf(stderr,"-------------------- begin P_Enum -------------------\n");
01236   fprintf(stderr,"Calling P_Enum with pos = %d\n",pos);
01237 #endif
01238   
01239   /* Xavier Redon modification: handle the case when there is no parameter */
01240   if(nb_param==0) {
01241     res=new_enode(polynomial,1,0);
01242     value_set_si(res->arr[0].d,1);
01243     value_init(res->arr[0].x.n);
01244     count_points(1,L,context,&res->arr[0].x.n);
01245     return res;
01246   }
01247   
01248   /* Initialize all the 'Value' variables */
01249   value_init(n); value_init(g); value_init(nLB);
01250   value_init(nUB); value_init(nlcm); value_init(noff);
01251   value_init(nexp); value_init(k1); value_init(nm);
01252   value_init(hdv); value_init(k); value_init(tmp);
01253   value_init(lcm_copy);
01254 
01255   if( value_zero_p(lcm[pos-1]) )
01256   {
01257                 hdim = 1;
01258                 value_set_si( lcm_copy, 1 );
01259   }
01260   else
01261   {
01262           /* hdim is the degree of the polynomial + 1 */
01263                 hdim = dim-nb_param+1;          /* homogenous dim w/o parameters */
01264                 value_assign( lcm_copy, lcm[pos-1] );
01265   }  
01266 
01267   /* code to limit generation of equations to valid parameters only */
01268   /*----------------------------------------------------------------*/
01269   flag = lower_upper_bounds(pos,LQ,&context[dim-nb_param],&nLB,&nUB);
01270   if (flag & LB_INFINITY) {
01271     if (!(flag & UB_INFINITY)) {
01272       
01273       /* Only an upper limit: set lower limit */
01274       /* Compute nLB such that (nUB-nLB+1) >= (hdim*lcm) */
01275       value_sub_int(nLB,nUB,1);
01276       value_set_si(hdv,hdim);
01277       value_multiply(tmp,hdv,lcm_copy);
01278       value_substract(nLB,nLB,tmp);
01279       if(value_pos_p(nLB))
01280         value_set_si(nLB,0);
01281     }
01282     else {
01283       value_set_si(nLB,0);
01284       
01285       /* No upper nor lower limit: set lower limit to 0 */
01286       value_set_si(hdv,hdim);
01287       value_multiply(nUB,hdv,lcm_copy);
01288       value_add_int(nUB,nUB,1);
01289     }
01290   }
01291 
01292   /* if (nUB-nLB+1) < (hdim*lcm) then we have more unknowns than equations */
01293   /* We can: 1. Find more equations by changing the context parameters, or */
01294   /* 2. Assign extra unknowns values in such a way as to simplify result.  */
01295   /* Think about ways to scan parameter space to get as much info out of it*/
01296   /* as possible.                                                          */
01297   
01298 #ifdef REDUCE_DEGREE
01299   if (pos==1 && (flag & UB_INFINITY)==0) {
01300     /* only for finite upper bound on first parameter */
01301     /* NOTE: only first parameter because subsequent parameters may
01302        be artificially limited by the choice of the first parameter */
01303    
01304 #ifdef EDEBUG
01305     fprintf(stderr,"*************** n **********\n");
01306     value_print(stderr,VALUE_FMT,n);
01307     fprintf(stderr,"\n");
01308 #endif
01309 
01310     value_substract(n,nUB,nLB);
01311     value_increment(n,n);
01312     
01313 #ifdef EDEBUG 
01314     value_print(stderr,VALUE_FMT,n);
01315     fprintf(stderr,"\n*************** n ************\n");
01316 #endif 
01317 
01318     /* Total number of samples>0 */
01319     if(value_neg_p(n))
01320       i=0;
01321     else {
01322       value_modulus(tmp,n,lcm_copy);
01323       if(value_notzero_p(tmp)) {
01324                   value_division(tmp,n,lcm_copy);
01325                   value_increment(tmp,tmp);
01326                   i = VALUE_TO_INT(tmp);
01327       }
01328       else {
01329                   value_division(tmp,n,lcm_copy);
01330                   i =   VALUE_TO_INT(tmp);      /* ceiling of n/lcm */
01331       }
01332     }
01333 
01334 #ifdef EDEBUG 
01335     value_print(stderr,VALUE_FMT,n);
01336     fprintf(stderr,"\n*************** n ************\n");
01337 #endif 
01338     
01339     /* Reduce degree of polynomial based on number of sample points */
01340     if (i < hdim){
01341       hdim=i;
01342       
01343 #ifdef EDEBUG4
01344       fprintf(stdout,"Parameter #%d: LB=",pos);
01345       value_print(stdout,VALUE_FMT,nLB);
01346       fprintf(stdout," UB=");
01347       value_print(stdout,VALUE_FMT,nUB);
01348       fprintf(stdout," lcm=");
01349       value_print(stdout,VALUE_FMT,lcm_copy);
01350       fprintf(stdout," degree reduced to %d\n",hdim-1);
01351 #endif
01352       
01353     }
01354   }
01355 #endif /* REDUCE_DEGREE */
01356   
01357   /* hdim is now set */  
01358   /* allocate result structure */
01359   res = new_enode(polynomial,hdim,pos);
01360   for (i=0;i<hdim;i++) {
01361     int l;
01362     l = VALUE_TO_INT(lcm_copy);
01363     res->arr[i].x.p = new_enode(periodic,l,pos);
01364   }
01365   
01366   /* Utility arrays */
01367   A = Matrix_Alloc(hdim, 2*hdim+1);                     /* used for Gauss */
01368   B = new_enode(evector, hdim, 0);
01369   C = new_enode(evector, hdim, 0);
01370   
01371   /*----------------------------------------------------------------*/
01372   /*                                                                */
01373   /*                                0<-----+---k--------->          */
01374   /* |---------noff----------------->-nlcm->-------lcm---->         */
01375   /* |--- . . . -----|--------------|------+-------|------+-------|-*/
01376   /* 0          (q-1)*lcm         q*lcm    |  (q+1)*lcm   |         */
01377   /*                                      nLB          nLB+lcm      */
01378   /*                                                                */
01379   /*----------------------------------------------------------------*/
01380   if(value_neg_p(nLB)) {
01381     value_modulus(nlcm,nLB,lcm_copy);
01382     value_addto(nlcm,nlcm,lcm_copy);
01383   }
01384   else {
01385     value_modulus(nlcm,nLB,lcm_copy);
01386   }
01387   
01388   /* noff is a multiple of lcm */
01389   value_substract(noff,nLB,nlcm);
01390   value_addto(tmp,lcm_copy,nlcm);
01391   for (value_assign(k,nlcm);value_lt(k,tmp);value_increment(k,k)) {
01392  
01393 #ifdef EDEBUG
01394     fprintf(stderr,"Finding ");
01395     value_print(stderr,VALUE_FMT,k);
01396     fprintf(stderr,"-th elements of periodic coefficients\n");
01397 #endif
01398 
01399     value_set_si(hdv,hdim);
01400     value_multiply(nm,hdv,lcm_copy);
01401     value_addto(nm,nm,nLB);
01402     i=0;
01403     for (value_addto(n,k,noff); value_lt(n,nm); value_addto(n,n,lcm_copy),i++) {
01404       
01405       /* n == i*lcm + k + noff;   */
01406       /* nlcm <= k < nlcm+lcm     */
01407       /* n mod lcm == k1 */
01408       
01409 #ifdef ALL_OVERFLOW_WARNINGS
01410       if (((flag & UB_INFINITY)==0) && value_gt(n,nUB)) {
01411         fprintf(stdout,"Domain Overflow: Parameter #%d:",pos);
01412         fprintf(stdout,"nLB=");
01413         value_print(stdout,VALUE_FMT,nLB);
01414         fprintf(stdout," n=");
01415         value_print(stdout,VALUE_FMT,n);
01416         fprintf(stdout," nUB=");
01417         value_print(stdout,VALUE_FMT,nUB);
01418         fprintf(stdout,"\n");
01419       }
01420 #else
01421       if (overflow_warning_flag && ((flag & UB_INFINITY)==0) &&
01422           value_gt(n,nUB)) {
01423         fprintf(stdout,"\nWARNING: Parameter Domain Overflow.");
01424         fprintf(stdout," Result may be incorrect on this domain.\n");
01425         overflow_warning_flag = 0;
01426       }
01427 #endif
01428       
01429       /* Set parameter to n */
01430       value_assign(context[dim-nb_param+pos],n);
01431       
01432 #ifdef EDEBUG1
01433         if( param_name )
01434         {
01435       fprintf(stderr,"%s = ",param_name[pos-1]);
01436       value_print(stderr,VALUE_FMT,n);
01437       fprintf(stderr," (hdim=%d, lcm[%d]=",hdim,pos-1);
01438       value_print(stderr,VALUE_FMT,lcm_copy);
01439       fprintf(stderr,")\n");
01440         }
01441         else
01442         {
01443       fprintf(stderr,"P%d = ",pos);
01444       value_print(stderr,VALUE_FMT,n);
01445       fprintf(stderr," (hdim=%d, lcm[%d]=",hdim,pos-1);
01446       value_print(stderr,VALUE_FMT,lcm_copy);
01447       fprintf(stderr,")\n");      
01448         }
01449 #endif
01450       
01451       /* Setup B vector */
01452       if (pos==nb_param) {
01453         
01454 #ifdef EDEBUG   
01455         fprintf(stderr,"Yes\n");
01456 #endif 
01457         
01458         /* call count */
01459         /* count can only be called when the context is fully specified */
01460         value_set_si(B->arr[i].d,1);
01461         value_init(B->arr[i].x.n);
01462         count_points(1,L,context,&B->arr[i].x.n);
01463         
01464 #ifdef EDEBUG3
01465         for (j=1; j<pos; j++) fputs("   ",stdout);
01466         fprintf(stdout, "E(");
01467         for (j=1; j<nb_param; j++) {
01468           value_print(stdout,VALUE_FMT,context[dim-nb_param+j]);
01469           fprintf(stdout,",");
01470         }  
01471         value_print(stdout,VALUE_FMT,n);
01472         fprintf(stdout,") = ");
01473         value_print(stdout,VALUE_FMT,B->arr[i].x.n);
01474         fprintf(stdout," =");
01475 #endif
01476         
01477       }
01478       else {    /* count is a function of other parameters */
01479         /* call P_Enum recursively */
01480         value_set_si(B->arr[i].d,0);
01481         B->arr[i].x.p = P_Enum(L,LQ->next,context,pos+1,nb_param,dim,lcm,param_name);
01482         
01483 #ifdef EDEBUG3
01484         if( param_name )
01485         {
01486                 for (j=1; j<pos; j++)
01487                   fputs("   ", stdout);
01488                 fprintf(stdout, "E(");
01489                 for (j=1; j<=pos; j++) {
01490                   value_print(stdout,VALUE_FMT,context[dim-nb_param+j]);
01491                   fprintf(stdout,",");
01492                 }         
01493                 for (j=pos+1; j<nb_param; j++)
01494                   fprintf(stdout,"%s,",param_name[j]);
01495                 fprintf(stdout,"%s) = ",param_name[j]);
01496                 print_enode(stdout,B->arr[i].x.p,param_name);
01497                 fprintf(stdout," =");
01498         }
01499 #endif
01500         
01501       }
01502       
01503       /* Create i-th equation*/
01504       /* K_0 + K_1 * n**1 + ... + K_dim * n**dim | Identity Matrix */
01505       
01506       value_set_si(A->p[i][0],0);  /* status bit = equality */
01507       value_set_si(nexp,1);
01508       for (j=1;j<=hdim;j++) {
01509         value_assign(A->p[i][j],nexp);
01510         value_set_si(A->p[i][j+hdim],0);
01511         
01512 #ifdef EDEBUG3
01513         fprintf(stdout," + ");
01514         value_print(stdout,VALUE_FMT,nexp);
01515         fprintf(stdout," c%d",j);
01516 #endif  
01517         value_multiply(nexp,nexp,n);
01518       }
01519       
01520 #ifdef EDEBUG3
01521       fprintf(stdout, "\n");
01522 #endif
01523       
01524       value_set_si(A->p[i][i+1+hdim],1);      
01525     }
01526     
01527     /* Assertion check */
01528     if (i!=hdim)
01529       fprintf(stderr, "P_Enum: ?expecting i==hdim\n");
01530     
01531 #ifdef EDEBUG
01532         if( param_name )
01533         {
01534     fprintf(stderr,"B (enode) =\n");
01535     print_enode(stderr,B,param_name);
01536         }
01537     fprintf(stderr,"A (Before Gauss) =\n");
01538     Matrix_Print(stderr,P_VALUE_FMT,A);
01539 #endif
01540     
01541     /* Solve hdim (=dim+1) equations with hdim unknowns, result in CNT */
01542     rank = Gauss(A,hdim,2*hdim);
01543     
01544 #ifdef EDEBUG
01545     fprintf(stderr,"A (After Gauss) =\n");
01546     Matrix_Print(stderr,P_VALUE_FMT,A);
01547 #endif
01548     
01549     /* Assertion check */
01550     if (rank!=hdim) {
01551       fprintf(stderr, "P_Enum: ?expecting rank==hdim\n");
01552     }
01553     
01554     /* if (rank < hdim) then set the last hdim-rank coefficients to ? */
01555     /* if (rank == 0) then set all coefficients to 0 */    
01556     /* copy result as k1-th element of periodic numbers */
01557     if(value_lt(k,lcm_copy))
01558       value_assign(k1,k);
01559     else
01560       value_substract(k1,k,lcm_copy);
01561     
01562     for (i=0; i<rank; i++) {
01563       
01564       /* Set up coefficient vector C from i-th row of inverted matrix */
01565       for (j=0; j<rank; j++) {
01566         Gcd(A->p[i][i+1],A->p[i][j+1+hdim],&g);
01567         value_division(C->arr[j].d,A->p[i][i+1],g);
01568         value_init(C->arr[j].x.n);
01569         value_division(C->arr[j].x.n,A->p[i][j+1+hdim],g);
01570       }
01571       
01572 #ifdef EDEBUG
01573         if( param_name )
01574         {
01575       fprintf(stderr, "C (enode) =\n");
01576       print_enode(stderr, C, param_name);
01577         }
01578 #endif
01579       
01580       /* The i-th enode is the lcm-periodic coefficient of term n**i */
01581       edot(B,C,&(res->arr[i].x.p->arr[VALUE_TO_INT(k1)]));
01582       
01583 #ifdef EDEBUG
01584         if( param_name )
01585         {
01586       fprintf(stderr, "B.C (evalue)=\n");
01587       print_evalue(stderr,&(res->arr[i].x.p->arr[VALUE_TO_INT(k1)]),param_name );
01588       fprintf(stderr,"\n");
01589         }
01590 #endif
01591       
01592     }
01593     value_addto(tmp,lcm_copy,nlcm);
01594   }
01595   
01596 #ifdef EDEBUG
01597         if( param_name )
01598         {
01599           fprintf(stderr,"res (enode) =\n");
01600           print_enode(stderr,res,param_name);
01601         }
01602   fprintf(stderr, "-------------------- end P_Enum -----------------------\n");
01603 #endif
01604   
01605   /* Reset context */
01606   value_set_si(context[dim-nb_param+pos],0);
01607   
01608   /* Release memory */
01609   Matrix_Free(A);
01610   free(B);
01611   free(C);
01612   
01613   /* Clear all the 'Value' variables */
01614   value_clear(n); value_clear(g); value_clear(nLB);
01615   value_clear(nUB); value_clear(nlcm); value_clear(noff);
01616   value_clear(nexp); value_clear(k1); value_clear(nm);
01617   value_clear(hdv); value_clear(k); value_clear(tmp);
01618   value_clear(lcm_copy);
01619   return res;
01620 } /* P_Enum */
01621 
01622 /*----------------------------------------------------------------*/
01623 /* Scan_Vertices(PP, Q, CT)                                       */
01624 /*    PP : ParamPolyhedron                                        */
01625 /*    Q : Domain                                                  */
01626 /*    CT : Context transformation matrix                          */
01627 /*    lcm : lcm array (output)                                    */
01628 /*    nbp : number of parameters                                  */
01629 /*    param_name : name of the parameters                         */
01630 /*----------------------------------------------------------------*/
01631 static void Scan_Vertices(Param_Polyhedron *PP,Param_Domain *Q,Matrix *CT,
01632    Value *lcm, int nbp, char **param_name )
01633 {
01634   Param_Vertices *V;
01635   int i, j, ix, l, np;
01636   unsigned bx;
01637   Value k,m1;
01638 
01639   /* Compute the denominator of P */
01640   /* lcm = Least Common Multiple of the denominators of the vertices of P */
01641   /* and print the vertices */  
01642 
01643   value_init(k); value_init(m1);
01644   for( np=0 ; np<nbp ; np++ )
01645     value_set_si( lcm[np], 0 );
01646 
01647         if( param_name )
01648           fprintf(stdout,"Vertices:\n");
01649 
01650   for(i=0,ix=0,bx=MSB,V=PP->V; V && i<PP->nbV; i++,V=V->next) {
01651     if (Q->F[ix] & bx) {
01652                 if( param_name )
01653                 {
01654         if(CT) {
01655                                 Matrix *v;
01656                                 v = VertexCT(V->Vertex,CT);
01657                                 Print_Vertex(stdout,v,param_name);
01658                                 Matrix_Free(v);
01659         }
01660         else
01661                                 Print_Vertex(stdout,V->Vertex,param_name);
01662                 fprintf(stdout,"\n");
01663                 }
01664 
01665       for(j=0;j<V->Vertex->NbRows;j++) {
01666                         /* A matrix */
01667                         for( l=0 ; l<V->Vertex->NbColumns-2 ; l++ )
01668                         {
01669                                 if( value_notzero_p(V->Vertex->p[j][l]) )
01670                                 {
01671                                         Gcd(V->Vertex->p[j][V->Vertex->NbColumns-1],V->Vertex->p[j][l], &m1);
01672                                         value_division(k,V->Vertex->p[j][V->Vertex->NbColumns-1],m1);
01673                                         if( value_notzero_p(lcm[l]) )
01674                                         {
01675                                                 /* lcm[l] = lcm[l] * k / gcd(k,lcm[l]) */
01676                                                 if (value_notzero_p(k) && value_notone_p(k))
01677                                                 {
01678                                                   Gcd(lcm[l],k,&m1);
01679                                                   value_division(k,k,m1);
01680                                                   value_multiply(lcm[l],lcm[l],k);
01681                                                 }
01682                                         }
01683                                         else
01684                                         {
01685                                                 value_assign(lcm[l],k);
01686                                         }
01687                                 }
01688                         }
01689       }
01690     }
01691     NEXT(ix,bx);
01692   }
01693   value_clear(k); value_clear(m1);
01694 } /* Scan_Vertices */
01695 
01696 /**
01697 
01698 Procedure to count points in a non-parameterized polytope.
01699 
01700 @param P       Polyhedron to count
01701 @param C       Parameter Context domain
01702 @param CT      Matrix to transform context to original
01703 @parma CEq     additionnal equalities in context
01704 @param MAXRAYS workspace size
01705 @param param_name parameter names
01706 
01707 */
01708 Enumeration *Enumerate_NoParameters(Polyhedron *P,Polyhedron *C,Matrix *CT,Polyhedron *CEq,unsigned MAXRAYS,char **param_name) {
01709   
01710     Polyhedron *L;
01711     Enumeration *res;
01712     Value *context,tmp;
01713     int j;
01714     int hdim = P->Dimension + 1;
01715     int r,i;
01716   
01717     /* Create a context vector size dim+2 */
01718     context = (Value *) malloc((hdim+1)*sizeof(Value));
01719     for (j=0;j<= hdim;j++) 
01720         value_init(context[j]);
01721     value_init(tmp);
01722   
01723     res = (Enumeration *)malloc(sizeof(Enumeration));
01724     res->next = NULL;
01725     res->ValidityDomain = Universe_Polyhedron(0);  /* no parameters */
01726     value_init(res->EP.d);
01727     value_set_si(res->EP.d,0);  
01728     L = Polyhedron_Scan(P,res->ValidityDomain,MAXRAYS);
01729   
01730 #ifdef EDEBUG2
01731     fprintf(stderr, "L = \n");
01732     Polyhedron_Print(stderr, P_VALUE_FMT, L);
01733 #endif
01734   
01735     if(CT) {
01736         Polyhedron *Dt;
01737     
01738         /* Add parameters to validity domain */
01739         Dt = Polyhedron_Preimage(res->ValidityDomain,CT,MAXRAYS);
01740         Polyhedron_Free(res->ValidityDomain);
01741         res->ValidityDomain = DomainIntersection(Dt,CEq,MAXRAYS);
01742         Polyhedron_Free(Dt);
01743     }
01744   
01745         if( param_name )
01746         {
01747     fprintf(stdout,"---------------------------------------\n");
01748     fprintf(stdout,"Domain:\n");
01749     Print_Domain(stdout,res->ValidityDomain, param_name);
01750   
01751     /* Print the vertices */
01752     printf("Vertices:\n");
01753     for(r=0;r<P->NbRays;++r) {
01754         if(value_zero_p(P->Ray[r][0]))
01755             printf("(line) ");
01756         printf("[");
01757         if (P->Dimension > 0)
01758             value_print(stdout,P_VALUE_FMT,P->Ray[r][1]);
01759         for(i=1;i<P->Dimension;i++) {
01760             printf(", ");
01761             value_print(stdout,P_VALUE_FMT,P->Ray[r][i+1]);
01762         }  
01763         printf("]");
01764         if(value_notone_p(P->Ray[r][P->Dimension+1])) {
01765             printf("/");
01766             value_print(stdout,P_VALUE_FMT, P->Ray[r][P->Dimension+1]);
01767         }
01768         printf("\n");
01769     }
01770         }
01771   if (emptyQ(P)) {
01772     res->EP.x.p = new_enode(polynomial,1,0);;
01773     value_set_si(res->EP.x.p->arr[0].d, 1);
01774     value_init(res->EP.x.p->arr[0].x.n);
01775     value_set_si(res->EP.x.p->arr[0].x.n, 0);
01776   } else if (!L) {
01777     /* Non-empty zero-dimensional domain */
01778     res->EP.x.p = new_enode(polynomial,1,0);;
01779     value_set_si(res->EP.x.p->arr[0].d, 1);
01780     value_init(res->EP.x.p->arr[0].x.n);
01781     value_set_si(res->EP.x.p->arr[0].x.n, 1);
01782   } else {
01783     CATCH(overflow_error) {
01784         fprintf(stderr,"Enumerate: arithmetic overflow error.\n");
01785         fprintf(stderr,"You should rebuild PolyLib using GNU-MP or increasing the size of integers.\n");
01786         overflow_warning_flag = 0;
01787         assert(overflow_warning_flag);
01788     }
01789     TRY {
01790     
01791         Vector_Set(context,0,(hdim+1));
01792     
01793         /* Set context[hdim] = 1  (the constant) */
01794         value_set_si(context[hdim],1);
01795         value_set_si(tmp,1);
01796         res->EP.x.p = P_Enum(L,NULL,context,1,0,hdim-1,&tmp,param_name);
01797         UNCATCH(overflow_error);
01798     }
01799   }
01800   
01801     Domain_Free(L);
01802   
01803     /*  **USELESS, there are no references to parameters in res**
01804         if( CT )
01805         addeliminatedparams_evalue(&res->EP, CT);
01806     */
01807   
01808         if( param_name )
01809         {
01810     fprintf(stdout,"\nEhrhart Polynomial:\n");
01811     print_evalue(stdout,&res->EP,param_name);
01812     fprintf(stdout, "\n");
01813         }
01814 
01815     value_clear(tmp);
01816     for (j=0;j<= hdim;j++) 
01817         value_clear(context[j]);  
01818     free(context);
01819     return(res);
01820 } /* Enumerate_NoParameters */
01821 
01822 /**  Procedure to count points in a parameterized polytope.
01823 
01824 @param P Polyhedron to enumerate
01825 @param C Context Domain
01826 @param MAXRAYS size of workspace
01827 @param param_name parameter names (array of strings), may be NULL
01828 @return a list of validity domains + evalues EP
01829 
01830 */
01831 Enumeration *Polyhedron_Enumerate(Polyhedron *Pi,Polyhedron *C,unsigned MAXRAYS,char **param_name)
01832 {
01833   Polyhedron *L, *CQ, *CQ2, *LQ, *U, *CEq, *rVD, *P;
01834   Matrix *CT;
01835   Param_Polyhedron *PP;
01836   Param_Domain   *Q;
01837   int i,hdim, dim, nb_param, np;
01838   Value *lcm, *m1, hdv;
01839   Value *context;
01840   Enumeration *en, *res;
01841 
01842   res = NULL;
01843   P = Pi;
01844 
01845 #ifdef EDEBUG2
01846   fprintf(stderr,"C = \n");
01847   Polyhedron_Print(stderr,P_VALUE_FMT,C);
01848   fprintf(stderr,"P = \n");
01849   Polyhedron_Print(stderr,P_VALUE_FMT,P);
01850 #endif
01851 
01852   hdim          = P->Dimension + 1;
01853   dim           = P->Dimension;
01854   nb_param      = C->Dimension;
01855 
01856   /* Don't call Polyhedron2Param_Domain if there are no parameters */
01857   if(nb_param == 0) {
01858     
01859     return(Enumerate_NoParameters(P,C,NULL,NULL,MAXRAYS,param_name));  
01860   }
01861   if(nb_param == dim) {
01862     res = (Enumeration *)malloc(sizeof(Enumeration));
01863     res->next = 0;
01864     res->ValidityDomain = DomainIntersection(P,C,MAXRAYS);
01865     value_init(res->EP.d);
01866     value_init(res->EP.x.n);
01867     value_set_si(res->EP.d,1);
01868     value_set_si(res->EP.x.n,1);
01869     if( param_name ) {
01870       fprintf(stdout,"---------------------------------------\n");
01871       fprintf(stdout,"Domain:\n");
01872       Print_Domain(stdout,res->ValidityDomain, param_name);
01873       fprintf(stdout,"\nEhrhart Polynomial:\n");
01874       print_evalue(stdout,&res->EP,param_name);
01875       fprintf(stdout, "\n");
01876     }
01877     return res;
01878   }
01879   PP = Polyhedron2Param_SimplifiedDomain(&P,C,MAXRAYS,&CEq,&CT);
01880   if(!PP) {
01881         if( param_name )
01882     fprintf(stdout, "\nEhrhart Polynomial:\nNULL\n");
01883 
01884     return(NULL);
01885   }
01886   
01887   /* CT : transformation matrix to eliminate useless ("false") parameters */
01888   if(CT) {
01889     nb_param -= CT->NbColumns-CT->NbRows;
01890     dim  -= CT->NbColumns-CT->NbRows;
01891     hdim -= CT->NbColumns-CT->NbRows;
01892     
01893     /* Don't call Polyhedron2Param_Domain if there are no parameters */
01894     if(nb_param == 0)
01895     {
01896             res = Enumerate_NoParameters(P,C,CT,CEq,MAXRAYS,param_name);
01897             if( P != Pi )
01898                     Polyhedron_Free( P );
01899             return( res );
01900     }  
01901   }
01902 
01903   /* get memory for Values */
01904   lcm = (Value *)malloc( nb_param * sizeof(Value));
01905   m1  = (Value *)malloc( nb_param * sizeof(Value));
01906   /* Initialize all the 'Value' variables */
01907   for( np=0 ; np<nb_param ; np++ )
01908   {
01909     value_init(lcm[np]); value_init(m1[np]);
01910   }
01911   value_init(hdv);
01912 
01913   for(Q=PP->D;Q;Q=Q->next) {
01914     if(CT) {
01915       Polyhedron *Dt;
01916       CQ = Q->Domain;      
01917       Dt = Polyhedron_Preimage(Q->Domain,CT,MAXRAYS);
01918       rVD = DomainIntersection(Dt,CEq,MAXRAYS);
01919       
01920       /* if rVD is empty or too small in geometric dimension */
01921       if(!rVD || emptyQ(rVD) ||
01922          (rVD->Dimension-rVD->NbEq < Dt->Dimension-Dt->NbEq-CEq->NbEq)) {
01923         if(rVD)
01924           Polyhedron_Free(rVD);
01925         Polyhedron_Free(Dt);
01926         continue;               /* empty validity domain */
01927       }
01928       Polyhedron_Free(Dt);
01929     }
01930     else
01931       rVD = CQ = Q->Domain;    
01932     en = (Enumeration *)malloc(sizeof(Enumeration));
01933     en->next = res;
01934     res = en;
01935     res->ValidityDomain = rVD;
01936     
01937         if( param_name )
01938    {
01939          fprintf(stdout,"---------------------------------------\n");
01940     fprintf(stdout,"Domain:\n");
01941     
01942 #ifdef EPRINT_ALL_VALIDITY_CONSTRAINTS
01943     Print_Domain(stdout,res->ValidityDomain,param_name);
01944 #else    
01945     {
01946       Polyhedron *VD;
01947       VD = DomainSimplify(res->ValidityDomain,C,MAXRAYS);
01948       Print_Domain(stdout,VD,param_name);
01949       Domain_Free(VD);
01950     }
01951 #endif /* EPRINT_ALL_VALIDITY_CONSTRAINTS */
01952         }
01953     
01954     overflow_warning_flag = 1;
01955     
01956     /* Scan the vertices and compute lcm */
01957     Scan_Vertices(PP,Q,CT,lcm,nb_param,param_name);
01958     
01959 #ifdef EDEBUG2
01960     fprintf(stderr,"Denominator = ");
01961     for( np=0;np<nb_param;np++)
01962       value_print(stderr,P_VALUE_FMT,lcm[np]);
01963     fprintf(stderr," and hdim == %d \n",hdim);
01964 #endif
01965     
01966 #ifdef EDEBUG2
01967     fprintf(stderr,"CQ = \n");
01968     Polyhedron_Print(stderr,P_VALUE_FMT,CQ);
01969 #endif
01970 
01971     /* Before scanning, add constraints to ensure at least hdim*lcm */
01972     /* points in every dimension */
01973     value_set_si(hdv,hdim-nb_param);
01974 
01975     for( np=0;np<nb_param;np++)
01976     {
01977                 if( value_notzero_p(lcm[np]) )
01978                         value_multiply(m1[np],hdv,lcm[np]);
01979                 else
01980                         value_set_si(m1[np],1);
01981     }
01982 
01983 #ifdef EDEBUG2 
01984     fprintf(stderr,"m1 == ");
01985     for( np=0;np<nb_param;np++)
01986       value_print(stderr,P_VALUE_FMT,m1[np]);
01987     fprintf(stderr,"\n");
01988 #endif 
01989 
01990     CATCH(overflow_error) {
01991       fprintf(stderr,"Enumerate: arithmetic overflow error.\n");
01992       CQ2 = NULL;
01993     }
01994     TRY {      
01995       CQ2 = Polyhedron_Preprocess(CQ,m1,MAXRAYS);
01996 
01997 #ifdef EDEBUG2
01998       fprintf(stderr,"After preprocess, CQ2 = ");
01999       Polyhedron_Print(stderr,P_VALUE_FMT,CQ2);
02000 #endif
02001       
02002       UNCATCH(overflow_error);
02003     }
02004     
02005     /* Vin100, Feb 2001 */
02006     /* in case of degenerate, try to find a domain _containing_ CQ */
02007     if ((!CQ2 || emptyQ(CQ2)) && CQ->NbBid==0) {
02008       int r;
02009       
02010 #ifdef EDEBUG2
02011       fprintf(stderr,"Trying to call Polyhedron_Preprocess2 : CQ = \n");
02012       Polyhedron_Print(stderr,P_VALUE_FMT,CQ);
02013 #endif
02014       
02015       /* Check if CQ is bounded */
02016       for(r=0;r<CQ->NbRays;r++) {
02017         if(value_zero_p(CQ->Ray[r][0]) ||
02018            value_zero_p(CQ->Ray[r][CQ->Dimension+1]))
02019           break;
02020       }
02021       if(r==CQ->NbRays) {
02022         
02023         /* ok, CQ is bounded */
02024         /* now find if CQ is contained in a hypercube of size m1 */
02025         CQ2 = Polyhedron_Preprocess2(CQ,m1,lcm,MAXRAYS);
02026       }
02027     }
02028     if (!CQ2 || emptyQ(CQ2)) {
02029 #ifdef EDEBUG2
02030       fprintf(stderr,"Degenerate.\n");
02031 #endif
02032       fprintf(stdout,"Degenerate Domain. Can not continue.\n");
02033       value_init(res->EP.d);
02034       value_init(res->EP.x.n);
02035       value_set_si(res->EP.d,1);
02036       value_set_si(res->EP.x.n,-1);
02037     }
02038     else {
02039       
02040 #ifdef EDEBUG2
02041       fprintf(stderr,"CQ2 = \n");
02042       Polyhedron_Print(stderr,P_VALUE_FMT,CQ2);
02043       if( ! PolyhedronIncludes(CQ, CQ2) )
02044         fprintf( stderr,"CQ does not include CQ2 !\n");
02045       else
02046         fprintf( stderr,"CQ includes CQ2.\n");
02047       if( ! PolyhedronIncludes(res->ValidityDomain, CQ2) )
02048         fprintf( stderr,"CQ2 is *not* included in validity domain !\n");
02049       else
02050         fprintf( stderr,"CQ2 is included in validity domain.\n");
02051 #endif
02052       
02053       /* L is used in counting the number of points in the base cases */
02054       L = Polyhedron_Scan(P,CQ,MAXRAYS);
02055       U = Universe_Polyhedron(0);
02056       
02057       /* LQ is used to scan the parameter space */
02058       LQ = Polyhedron_Scan(CQ2,U,MAXRAYS); /* bounds on parameters */
02059       Domain_Free(U);
02060       if(CT)    /* we did compute another Q->Domain */
02061         Domain_Free(CQ);
02062       
02063       /* Else, CQ was Q->Domain (used in res) */
02064       Domain_Free(CQ2);
02065       
02066 #ifdef EDEBUG2
02067       fprintf(stderr,"L = \n");
02068       Polyhedron_Print(stderr,P_VALUE_FMT,L);
02069       fprintf(stderr,"LQ = \n");
02070       Polyhedron_Print(stderr,P_VALUE_FMT,LQ);
02071 #endif
02072 #ifdef EDEBUG3
02073       fprintf(stdout,"\nSystem of Equations:\n");
02074 #endif
02075       
02076       value_init(res->EP.d);
02077       value_set_si(res->EP.d,0);
02078       
02079       /* Create a context vector size dim+2 */
02080       context = (Value *) malloc ((hdim+1)*sizeof(Value));  
02081       for(i=0;i<=(hdim);i++)
02082         value_init(context[i]);
02083       Vector_Set(context,0,(hdim+1));
02084       
02085       /* Set context[hdim] = 1  (the constant) */
02086       value_set_si(context[hdim],1);
02087       
02088       CATCH(overflow_error) {
02089                         fprintf(stderr,"Enumerate: arithmetic overflow error.\n");
02090                         fprintf(stderr,"You should rebuild PolyLib using GNU-MP or increasing the size of integers.\n");
02091                         overflow_warning_flag = 0;
02092                         assert(overflow_warning_flag);
02093         
02094       }
02095       TRY {
02096                         res->EP.x.p = P_Enum(L,LQ,context,1,nb_param,dim,lcm,param_name);
02097                         UNCATCH(overflow_error);        
02098       }
02099       
02100       for(i=0;i<=(hdim);i++)
02101         value_clear(context[i]);
02102       free(context);
02103       Domain_Free(L);
02104       Domain_Free(LQ);
02105       
02106 #ifdef EDEBUG5
02107         if( param_name )
02108    {
02109       fprintf(stdout,"\nEhrhart Polynomial (before simplification):\n");
02110       print_evalue(stdout,&res->EP,param_name);
02111         }
02112 #endif
02113 
02114       /* Try to simplify the result */
02115       reduce_evalue(&res->EP);
02116       
02117       /* Put back the original parameters into result */
02118       /* (equalities have been eliminated)            */
02119       if(CT) 
02120           addeliminatedparams_evalue(&res->EP,CT);
02121       
02122         if( param_name )
02123         {
02124       fprintf(stdout,"\nEhrhart Polynomial:\n");
02125       print_evalue(stdout,&res->EP, param_name);
02126       fprintf(stdout,"\n");
02127       /* sometimes the final \n lacks (when a single constant is printed) */
02128         }
02129       
02130     }
02131   }
02132 
02133   if( P != Pi )
02134         Polyhedron_Free( P );
02135   /* Clear all the 'Value' variables */
02136   for( np=0; np<nb_param ; np++ )
02137   {
02138     value_clear(lcm[np]); value_clear(m1[np]);
02139   }
02140   value_clear(hdv);
02141 
02142   return res;
02143 } /* Polyhedron_Enumerate */ 
02144 
02145 void Enumeration_Free(Enumeration *en)
02146 {
02147   Enumeration *ee;
02148 
02149   while( en )
02150   {
02151           free_evalue_refs( &(en->EP) );
02152           Polyhedron_Free( en->ValidityDomain );
02153           ee = en ->next;
02154           free( en );
02155           en = ee;
02156   }
02157 }
02158 
02159 
02160 

Generated on Mon Sep 12 14:48:28 2005 for polylib by doxygen 1.3.5