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

matrix.c

Go to the documentation of this file.
00001 /* matrix.c 
00002      COPYRIGHT
00003           Both this software and its documentation are
00004 
00005               Copyright 1993 by IRISA /Universite de Rennes I -
00006               France, Copyright 1995,1996 by BYU, Provo, Utah
00007                          all rights reserved.
00008 
00009           Permission is granted to copy, use, and distribute
00010           for any commercial or noncommercial purpose under the terms
00011           of the GNU General Public license, version 2, June 1991
00012           (see file : LICENSING).
00013 */
00014 
00015 #include <stdio.h>
00016 #include <stdlib.h>
00017 #include <string.h>
00018 #include <ctype.h>
00019 #include <polylib/polylib.h>
00020 
00021 #ifdef mac_os
00022   #define abs __abs
00023 #endif
00024 
00025 /* 
00026  * Allocate space for matrix dimensioned by 'NbRows X NbColumns'.
00027  */
00028 Matrix *Matrix_Alloc(unsigned NbRows,unsigned NbColumns) {
00029   
00030   Matrix *Mat;
00031   Value *p, **q;
00032   int i,j;
00033 
00034   Mat=(Matrix *)malloc(sizeof(Matrix));
00035   if(!Mat) {    
00036     errormsg1("Matrix_Alloc", "outofmem", "out of memory space");
00037     return 0;
00038   }
00039   Mat->NbRows=NbRows;
00040   Mat->NbColumns=NbColumns;
00041   if (NbRows==0 || NbColumns==0) {
00042       Mat->p = (Value **)0;
00043       Mat->p_Init= (Value *)0;
00044       Mat->p_Init_size = 0;
00045   } else {
00046       q = (Value **)malloc(NbRows * sizeof(*q));
00047       if(!q) {
00048         free(Mat);
00049         errormsg1("Matrix_Alloc", "outofmem", "out of memory space");
00050         return 0;
00051       }
00052       p = value_alloc(NbRows * NbColumns, &Mat->p_Init_size);
00053       if(!p) {
00054         free(q);
00055         free(Mat);
00056         errormsg1("Matrix_Alloc", "outofmem", "out of memory space");
00057         return 0;
00058       }
00059       Mat->p = q;
00060       Mat->p_Init = p;
00061       for (i=0;i<NbRows;i++) {
00062         *q++ = p;
00063         p += NbColumns;
00064       }
00065   }
00066   p = NULL;
00067   q = NULL;
00068 
00069   return Mat;
00070 } /* Matrix_Alloc */
00071 
00072 /* 
00073  * Free the memory space occupied by Matrix 'Mat' 
00074  */
00075 void Matrix_Free(Matrix *Mat)
00076 { 
00077   if (Mat->p_Init)
00078     value_free(Mat->p_Init, Mat->p_Init_size);
00079 
00080   if (Mat->p)
00081     free(Mat->p);
00082   free(Mat);
00083 
00084 } /* Matrix_Free */
00085 
00086 void Matrix_Extend(Matrix *Mat, unsigned NbRows)
00087 {
00088   Value *p, **q;
00089   int i,j;
00090 
00091   q = (Value **)realloc(Mat->p, NbRows * sizeof(*q));
00092   if(!q) {
00093     errormsg1("Matrix_Extend", "outofmem", "out of memory space");
00094     return;
00095   }
00096   Mat->p = q;
00097   if (Mat->p_Init_size < NbRows * Mat->NbColumns) {
00098     p = (Value *)realloc(Mat->p_Init, NbRows * Mat->NbColumns * sizeof(Value));
00099     if(!p) {
00100       errormsg1("Matrix_Extend", "outofmem", "out of memory space");
00101       return;
00102     }
00103     Mat->p_Init = p;
00104     Vector_Set(Mat->p_Init + Mat->NbRows*Mat->NbColumns, 0,
00105                Mat->p_Init_size - Mat->NbRows*Mat->NbColumns);
00106     for (i = Mat->p_Init_size; i < Mat->NbColumns*NbRows; ++i)
00107         value_init(Mat->p_Init[i]);
00108     Mat->p_Init_size = Mat->NbColumns*NbRows;
00109   } else
00110     Vector_Set(Mat->p_Init + Mat->NbRows*Mat->NbColumns, 0,
00111                (NbRows - Mat->NbRows) * Mat->NbColumns);
00112   for (i=0;i<NbRows;i++) {
00113     Mat->p[i] = Mat->p_Init + (i * Mat->NbColumns);
00114   }
00115   Mat->NbRows = NbRows;
00116 }
00117 
00118 /* 
00119  * Print the contents of the Matrix 'Mat'
00120  */
00121 void Matrix_Print(FILE *Dst, const char *Format, Matrix *Mat)
00122 {
00123   Value *p;
00124   int i, j;
00125   unsigned NbRows, NbColumns;
00126   
00127   fprintf(Dst,"%d %d\n", NbRows=Mat->NbRows, NbColumns=Mat->NbColumns);
00128   if (NbColumns ==0) {
00129     fprintf(Dst, "\n");
00130     return;
00131   }
00132   for (i=0;i<NbRows;i++) {
00133     p=*(Mat->p+i);
00134     for (j=0;j<NbColumns;j++) {
00135       if (!Format) {
00136         value_print(Dst," "P_VALUE_FMT" ",*p++);
00137       }
00138       else { 
00139         value_print(Dst,Format,*p++);
00140       } 
00141     }
00142     fprintf(Dst, "\n");
00143   }
00144 } /* Matrix_Print */
00145 
00146 /* 
00147  * Read the contents of the Matrix 'Mat' 
00148  */
00149 void Matrix_Read_Input(Matrix *Mat) {
00150   
00151   Value *p;
00152   int i,j,n;
00153   char *c, s[1024],str[1024];
00154   
00155   p = Mat->p_Init;
00156   for (i=0;i<Mat->NbRows;i++) {
00157     do {
00158       c = fgets(s, 1024, stdin);
00159       while(isspace(*c) && *c!='\n')
00160         ++c;
00161     } while(c && (*c =='#' || *c== '\n'));
00162     
00163     if (!c) {
00164       errormsg1( "Matrix_Read", "baddim", "not enough rows" );
00165       break;
00166     }
00167     for (j=0;j<Mat->NbColumns;j++) {
00168       if(!c || *c=='\n' || *c=='#') {
00169         errormsg1("Matrix_Read", "baddim", "not enough columns");
00170         break;
00171       }
00172       if (sscanf(c,"%s%n",str,&n) == 0) {
00173         errormsg1( "Matrix_Read", "baddim", "not enough columns" );
00174         break;
00175       }
00176       value_read(*(p++),str);
00177       c += n;
00178     }
00179   }
00180 } /* Matrix_Read_Input */
00181 
00182 /* 
00183  * Read the contents of the matrix 'Mat' from standard input. 
00184  * A '#' in the first column is a comment line 
00185  */
00186 Matrix *Matrix_Read(void) {
00187   
00188   Matrix *Mat;
00189   unsigned NbRows, NbColumns;
00190   char s[1024];
00191   
00192   if (fgets(s, 1024, stdin) == NULL)
00193     return NULL;
00194   while ((*s=='#' || *s=='\n') ||
00195          (sscanf(s, "%d %d", &NbRows, &NbColumns)<2)) {
00196     if (fgets(s, 1024, stdin) == NULL)
00197       return NULL;
00198   }
00199   Mat = Matrix_Alloc(NbRows,NbColumns);
00200   if(!Mat) {
00201     errormsg1("Matrix_Read", "outofmem", "out of memory space");
00202     return(NULL);
00203   }
00204   Matrix_Read_Input(Mat);
00205   return Mat;
00206 } /* Matrix_Read */
00207 
00208 /* 
00209  * Basic hermite engine 
00210  */
00211 static int hermite(Matrix *H,Matrix *U,Matrix *Q) {
00212   
00213   int nc, nr, i, j, k, rank, reduced, pivotrow;
00214   Value pivot,x,aux;
00215   Value *temp1, *temp2;
00216   
00217   /*                     T                     -1   T */
00218   /* Computes form: A = Q H  and U A = H  and U  = Q  */
00219   
00220   if (!H) { 
00221     errormsg1("Domlib", "nullH", "hermite: ? Null H");
00222     return -1;
00223   }
00224   nc = H->NbColumns;
00225   nr = H->NbRows;
00226   temp1 = (Value *) malloc(nc * sizeof(Value));
00227   temp2 = (Value *) malloc(nr * sizeof(Value));
00228   if (!temp1 ||!temp2) {
00229     errormsg1("Domlib", "outofmem", "out of memory space");
00230     return -1;
00231   }
00232   
00233   /* Initialize all the 'Value' variables */
00234   value_init(pivot); value_init(x); 
00235   value_init(aux);   
00236   for(i=0;i<nc;i++)
00237     value_init(temp1[i]);
00238   for(i=0;i<nr;i++)
00239     value_init(temp2[i]);
00240   
00241 #ifdef DEBUG
00242   fprintf(stderr,"Start  -----------\n");
00243   Matrix_Print(stderr,0,H);
00244 #endif
00245   for (k=0, rank=0; k<nc && rank<nr; k=k+1) {
00246     reduced = 1;        /* go through loop the first time */
00247 #ifdef DEBUG
00248     fprintf(stderr, "Working on col %d.  Rank=%d ----------\n", k+1, rank+1);
00249 #endif
00250     while (reduced) {
00251       reduced=0;
00252       
00253       /* 1. find pivot row */
00254       value_absolute(pivot,H->p[rank][k]);
00255       
00256       /* the kth-diagonal element */
00257       pivotrow = rank;
00258       
00259       /* find the row i>rank with smallest nonzero element in col k */
00260       for (i=rank+1; i<nr; i++) {
00261         value_absolute(x,H->p[i][k]);
00262         if (value_notzero_p(x) &&
00263             (value_lt(x,pivot) || value_zero_p(pivot))) {
00264           value_assign(pivot,x);
00265           pivotrow = i;
00266         }
00267       }
00268       
00269       /* 2. Bring pivot to diagonal (exchange rows pivotrow and rank) */
00270       if (pivotrow != rank) {
00271         Vector_Exchange(H->p[pivotrow],H->p[rank],nc);
00272         if (U)
00273           Vector_Exchange(U->p[pivotrow],U->p[rank],nr);
00274         if (Q)
00275           Vector_Exchange(Q->p[pivotrow],Q->p[rank],nr);
00276 
00277 #ifdef DEBUG
00278         fprintf(stderr,"Exchange rows %d and %d  -----------\n", rank+1, pivotrow+1);
00279         Matrix_Print(stderr,0,H);
00280 #endif
00281       }
00282       value_assign(pivot,H->p[rank][k]);        /* actual ( no abs() ) pivot */
00283       
00284       /* 3. Invert the row 'rank' if pivot is negative */
00285       if (value_neg_p(pivot)) {
00286         value_oppose(pivot,pivot); /* pivot = -pivot */
00287         for (j=0; j<nc; j++)
00288           value_oppose(H->p[rank][j],H->p[rank][j]);
00289         
00290         /* H->p[rank][j] = -(H->p[rank][j]); */
00291         if (U)
00292           for (j=0; j<nr; j++)
00293             value_oppose(U->p[rank][j],U->p[rank][j]);
00294         
00295         /* U->p[rank][j] = -(U->p[rank][j]); */
00296         if (Q)
00297           for (j=0; j<nr; j++)
00298             value_oppose(Q->p[rank][j],Q->p[rank][j]);
00299         
00300         /* Q->p[rank][j] = -(Q->p[rank][j]); */
00301 #ifdef DEBUG
00302         fprintf(stderr,"Negate row %d  -----------\n", rank+1);
00303         Matrix_Print(stderr,0,H);
00304 #endif
00305 
00306       }      
00307       if (value_notzero_p(pivot)) {
00308         
00309         /* 4. Reduce the column modulo the pivot */
00310         /*    This eventually zeros out everything below the */
00311         /*    diagonal and produces an upper triangular matrix */
00312         
00313         for (i=rank+1;i<nr;i++) {
00314           value_assign(x,H->p[i][k]);
00315           if (value_notzero_p(x)) {         
00316             value_modulus(aux,x,pivot);
00317             
00318             /* floor[integer division] (corrected for neg x) */
00319             if (value_neg_p(x) && value_notzero_p(aux)) {
00320               
00321               /* x=(x/pivot)-1; */
00322               value_division(x,x,pivot);
00323               value_decrement(x,x);
00324             }   
00325             else 
00326               value_division(x,x,pivot);
00327             for (j=0; j<nc; j++) {
00328               value_multiply(aux,x,H->p[rank][j]);
00329               value_subtract(H->p[i][j],H->p[i][j],aux);
00330             }
00331             
00332             /* U->p[i][j] -= (x * U->p[rank][j]); */
00333             if (U)
00334               for (j=0; j<nr; j++) {
00335                 value_multiply(aux,x,U->p[rank][j]);
00336                 value_subtract(U->p[i][j],U->p[i][j],aux);
00337               }
00338             
00339             /* Q->p[rank][j] += (x * Q->p[i][j]); */
00340             if (Q)
00341               for(j=0;j<nr;j++) {
00342                 value_addmul(Q->p[rank][j], x, Q->p[i][j]);
00343               }
00344             reduced = 1;
00345 
00346 #ifdef DEBUG
00347             fprintf(stderr,
00348                     "row %d = row %d - %d row %d -----------\n", i+1, i+1, x, rank+1);
00349             Matrix_Print(stderr,0,H);
00350 #endif
00351         
00352           } /* if (x) */
00353         } /* for (i) */
00354       } /* if (pivot != 0) */
00355     } /* while (reduced) */
00356     
00357     /* Last finish up this column */
00358     /* 5. Make pivot column positive (above pivot row) */
00359     /*    x should be zero for i>k */
00360     
00361     if (value_notzero_p(pivot)) {
00362       for (i=0; i<rank; i++) {
00363         value_assign(x,H->p[i][k]);
00364         if (value_notzero_p(x)) {         
00365           value_modulus(aux,x,pivot);
00366           
00367           /* floor[integer division] (corrected for neg x) */
00368           if (value_neg_p(x) && value_notzero_p(aux)) {
00369             value_division(x,x,pivot);
00370             value_decrement(x,x);
00371             
00372             /* x=(x/pivot)-1; */
00373           }
00374           else
00375             value_division(x,x,pivot);
00376           
00377           /* H->p[i][j] -= x * H->p[rank][j]; */
00378           for (j=0; j<nc; j++) {
00379             value_multiply(aux,x,H->p[rank][j]);
00380             value_subtract(H->p[i][j],H->p[i][j],aux);
00381           }
00382           
00383           /* U->p[i][j] -= x * U->p[rank][j]; */
00384           if (U)
00385             for (j=0; j<nr; j++) {
00386               value_multiply(aux,x,U->p[rank][j]);
00387               value_subtract(U->p[i][j],U->p[i][j],aux);
00388             }
00389           
00390           /* Q->p[rank][j] += x * Q->p[i][j]; */
00391           if (Q)
00392             for (j=0; j<nr; j++) {
00393               value_addmul(Q->p[rank][j], x, Q->p[i][j]);
00394             }  
00395 #ifdef DEBUG
00396           fprintf(stderr,
00397                   "row %d = row %d - %d row %d -----------\n", i+1, i+1, x, rank+1);
00398           Matrix_Print(stderr,0,H);
00399 #endif
00400         } /* if (x) */
00401       } /* for (i) */
00402       rank++;
00403     } /* if (pivot!=0) */
00404   } /* for (k) */
00405   
00406   /* Clear all the 'Value' variables */
00407   value_clear(pivot); value_clear(x); 
00408   value_clear(aux); 
00409   for(i=0;i<nc;i++)
00410     value_clear(temp1[i]);
00411   for(i=0;i<nr;i++)
00412     value_clear(temp2[i]);
00413   free(temp2);
00414   free(temp1);
00415   return rank;
00416 } /* Hermite */ 
00417 
00418 void right_hermite(Matrix *A,Matrix **Hp,Matrix **Up,Matrix **Qp) {
00419   
00420   Matrix *H, *Q, *U;
00421   int i, j, nr, nc, rank;
00422   Value tmp;
00423   
00424   /* Computes form: A = QH , UA = H */  
00425   nc = A->NbColumns;
00426   nr = A->NbRows;
00427   
00428   /* H = A */
00429   *Hp = H = Matrix_Alloc(nr,nc);
00430   if (!H) { 
00431     errormsg1("DomRightHermite", "outofmem", "out of memory space");
00432     return;
00433   }
00434   
00435   /* Initialize all the 'Value' variables */
00436   value_init(tmp);
00437   
00438   Vector_Copy(A->p_Init,H->p_Init,nr*nc);
00439   
00440   /* U = I */
00441   if (Up) {
00442     *Up = U = Matrix_Alloc(nr, nr);
00443     if (!U) {
00444       errormsg1("DomRightHermite", "outofmem", "out of memory space");
00445       value_clear(tmp);
00446       return;
00447     }
00448     Vector_Set(U->p_Init,0,nr*nr);             /* zero's */
00449     for(i=0;i<nr;i++)                          /* with diagonal of 1's */
00450       value_set_si(U->p[i][i],1);
00451   }
00452   else
00453     U = (Matrix *)0;
00454   
00455   /* Q = I */
00456   /* Actually I compute Q transpose... its easier */
00457   if (Qp) {
00458     *Qp = Q = Matrix_Alloc(nr,nr);
00459     if (!Q) {
00460       errormsg1("DomRightHermite", "outofmem", "out of memory space");
00461       value_clear(tmp);
00462       return;
00463     }
00464     Vector_Set(Q->p_Init,0,nr*nr);            /* zero's */
00465     for (i=0;i<nr;i++)                      /* with diagonal of 1's */
00466       value_set_si(Q->p[i][i],1);
00467   }
00468   else
00469     Q = (Matrix *)0;
00470   
00471   rank = hermite(H,U,Q);
00472   
00473   /* Q is returned transposed */ 
00474   /* Transpose Q */
00475   if (Q) {
00476     for (i=0; i<nr; i++) {
00477       for (j=i+1; j<nr; j++) {
00478         value_assign(tmp,Q->p[i][j]);
00479         value_assign(Q->p[i][j],Q->p[j][i] );
00480         value_assign(Q->p[j][i],tmp);
00481       }
00482     }
00483   }
00484   value_clear(tmp);
00485   return;
00486 } /* right_hermite */
00487 
00488 void left_hermite(Matrix *A,Matrix **Hp,Matrix **Qp,Matrix **Up) {
00489   
00490   Matrix *H, *HT, *Q, *U;
00491   int i, j, nc, nr, rank;
00492   Value tmp;
00493   
00494   /* Computes left form: A = HQ , AU = H , 
00495                         T    T T    T T   T
00496      using right form  A  = Q H  , U A = H */
00497   
00498   nr = A->NbRows;
00499   nc = A->NbColumns;
00500   
00501   /* HT = A transpose */
00502   HT = Matrix_Alloc(nc, nr);
00503   if (!HT) {
00504     errormsg1("DomLeftHermite", "outofmem", "out of memory space");
00505     return;
00506   }
00507   value_init(tmp);
00508   for (i=0; i<nr; i++)
00509     for (j=0; j<nc; j++)
00510       value_assign(HT->p[j][i],A->p[i][j]);
00511   
00512   /* U = I */
00513   if (Up) {
00514     *Up = U = Matrix_Alloc(nc,nc);
00515     if (!U) {
00516       errormsg1("DomLeftHermite", "outofmem", "out of memory space");
00517       value_clear(tmp);
00518       return;
00519     }
00520     Vector_Set(U->p_Init,0,nc*nc);            /* zero's */
00521     for (i=0;i<nc;i++)                        /* with diagonal of 1's */
00522       value_set_si(U->p[i][i],1);
00523   }
00524   else U=(Matrix *)0;
00525   
00526   /* Q = I */
00527   if (Qp) {
00528     *Qp = Q = Matrix_Alloc(nc, nc);
00529     if (!Q) {
00530       errormsg1("DomLeftHermite", "outofmem", "out of memory space");
00531       value_clear(tmp);
00532       return;
00533     }
00534     Vector_Set(Q->p_Init,0,nc*nc);            /* zero's */
00535     for (i=0;i<nc;i++)                        /* with diagonal of 1's */
00536       value_set_si(Q->p[i][i],1);
00537   }
00538   else Q=(Matrix *)0;
00539   rank = hermite(HT,U,Q);
00540   
00541   /* H = HT transpose */
00542   *Hp = H = Matrix_Alloc(nr,nc);
00543   if (!H) {
00544     errormsg1("DomLeftHermite", "outofmem", "out of memory space");
00545     value_clear(tmp);
00546     return;
00547   }
00548   for (i=0; i<nr; i++)
00549     for (j=0;j<nc;j++)
00550       value_assign(H->p[i][j],HT->p[j][i]);
00551   Matrix_Free(HT);
00552   
00553   /* Transpose U */
00554   if (U) {
00555     for (i=0; i<nc; i++) {
00556       for (j=i+1; j<nc; j++) {
00557         value_assign(tmp,U->p[i][j]);
00558         value_assign(U->p[i][j],U->p[j][i] );
00559         value_assign(U->p[j][i],tmp);
00560       }
00561     }
00562   }
00563   value_clear(tmp);
00564 } /* left_hermite */
00565 
00566 /*
00567  * Given a integer matrix 'Mat'(k x k), compute its inverse rational matrix 
00568  * 'MatInv' k x (k+1). The last column of each row in matrix MatInv is used 
00569  * to store the common denominator of the entries in a row. The output is 1,
00570  * if 'Mat' is non-singular (invertible), otherwise the output is 0. Note:: 
00571  * (1) Matrix 'Mat' is modified during the inverse operation.
00572  * (2) Matrix 'MatInv' must be preallocated before passing into this function.
00573  */
00574 int MatInverse(Matrix *Mat,Matrix *MatInv ) {
00575   
00576   int i, k, j, c;
00577   Value x, gcd, piv;
00578   Value m1,m2;
00579   
00580   if(Mat->NbRows != Mat->NbColumns) {
00581    fprintf(stderr,"Trying to invert a non-square matrix !\n");
00582     return 0;
00583   }
00584   
00585   /* Initialize all the 'Value' variables */
00586   value_init(x);  value_init(gcd); value_init(piv);
00587   value_init(m1); value_init(m2);
00588 
00589   k = Mat->NbRows; 
00590 
00591   /* Initialise MatInv */
00592   Vector_Set(MatInv->p[0],0,k*(k+1));
00593 
00594   /* Initialize 'MatInv' to Identity matrix form. Each diagonal entry is set*/
00595   /* to 1. Last column of each row (denominator of each entry in a row) is  */
00596   /* also set to 1.                                                         */ 
00597   for(i=0;i<k;++i) {
00598     value_set_si(MatInv->p[i][i],1);    
00599     value_set_si(MatInv->p[i][k],1);    /* denum */
00600   }  
00601   /* Apply Gauss-Jordan elimination method on the two matrices 'Mat' and  */
00602   /* 'MatInv' in parallel.                                                */
00603   for(i=0;i<k;++i) {
00604     
00605     /* Check if the diagonal entry (new pivot) is non-zero or not */
00606     if(value_zero_p(Mat->p[i][i])) {    
00607       
00608       /* Search for a non-zero pivot down the column(i) */
00609       for(j=i;j<k;++j)      
00610         if(value_notzero_p(Mat->p[j][i]))
00611           break;
00612       
00613       /* If no non-zero pivot is found, the matrix 'Mat' is non-invertible */
00614       /* Return 0.                                                         */
00615       if(j==k) {
00616         
00617         /* Clear all the 'Value' variables */
00618         value_clear(x);  value_clear(gcd); value_clear(piv);
00619         value_clear(m1); value_clear(m2);
00620         return 0;
00621       } 
00622       
00623       /* Exchange the rows, row(i) and row(j) so that the diagonal element */
00624       /* Mat->p[i][i] (pivot) is non-zero. Repeat the same operations on    */
00625       /* matrix 'MatInv'.                                                   */
00626       for(c=0;c<k;++c) {
00627 
00628         /* Interchange rows, row(i) and row(j) of matrix 'Mat'    */
00629         value_assign(x,Mat->p[j][c]);
00630         value_assign(Mat->p[j][c],Mat->p[i][c]);
00631         value_assign(Mat->p[i][c],x);
00632         
00633         /* Interchange rows, row(i) and row(j) of matrix 'MatInv' */
00634         value_assign(x,MatInv->p[j][c]);
00635         value_assign(MatInv->p[j][c],MatInv->p[i][c]);
00636         value_assign(MatInv->p[i][c],x);
00637       }
00638     }
00639     
00640     /* Make all the entries in column(i) of matrix 'Mat' zero except the */
00641     /* diagonal entry. Repeat the same sequence of operations on matrix  */
00642     /* 'MatInv'.                                                         */
00643     for(j=0;j<k;++j) {
00644       if (j==i) continue;                /* Skip the pivot */
00645       value_assign(x,Mat->p[j][i]);
00646       if(value_notzero_p(x)) {
00647         value_assign(piv,Mat->p[i][i]);
00648         value_gcd(gcd, x, piv);
00649         if (value_notone_p(gcd) ) {
00650           value_divexact(x, x, gcd);
00651           value_divexact(piv, piv, gcd);
00652         }
00653         for(c=((j>i)?i:0);c<k;++c) {
00654           value_multiply(m1,piv,Mat->p[j][c]);
00655           value_multiply(m2,x,Mat->p[i][c]);
00656           value_subtract(Mat->p[j][c],m1,m2); 
00657         }
00658         for(c=0;c<k;++c) {
00659           value_multiply(m1,piv,MatInv->p[j][c]);
00660           value_multiply(m2,x,MatInv->p[i][c]);
00661           value_subtract(MatInv->p[j][c],m1,m2);
00662         }
00663               
00664         /* Simplify row(j) of the two matrices 'Mat' and 'MatInv' by */
00665         /* dividing the rows with the common GCD.                     */
00666         Vector_Gcd(&MatInv->p[j][0],k,&m1);
00667         Vector_Gcd(&Mat->p[j][0],k,&m2);
00668         value_gcd(gcd, m1, m2);
00669         if(value_notone_p(gcd)) {
00670           for(c=0;c<k;++c) {
00671             value_divexact(Mat->p[j][c], Mat->p[j][c], gcd);
00672             value_divexact(MatInv->p[j][c], MatInv->p[j][c], gcd);
00673           }
00674         }
00675       }
00676     }
00677   }
00678   
00679   /* Simplify every row so that 'Mat' reduces to Identity matrix. Perform  */
00680   /* the same sequence of operations on the matrix 'MatInv'.               */
00681   for(j=0;j<k;++j) {
00682     value_assign(MatInv->p[j][k],Mat->p[j][j]);
00683     
00684     /* Make the last column (denominator of each entry) of every row greater */
00685     /* than zero.                                                            */
00686     Vector_Normalize_Positive(&MatInv->p[j][0],k+1,k); 
00687   }
00688   
00689   /* Clear all the 'Value' variables */
00690   value_clear(x);  value_clear(gcd); value_clear(piv);
00691   value_clear(m1); value_clear(m2);
00692 
00693   return 1;
00694 } /* Mat_Inverse */
00695 
00696 /*
00697  * Given (m x n) integer matrix 'X' and n x (k+1) rational matrix 'P', compute
00698  * the rational m x (k+1) rational matrix  'S'. The last column in each row of
00699  * the rational matrices is used to store the common denominator of elements
00700  * in a row.                              
00701  */
00702 void rat_prodmat(Matrix *S,Matrix *X,Matrix *P) {
00703   
00704   int i,j,k;
00705   int last_column_index = P->NbColumns - 1;
00706   Value lcm, old_lcm,gcd,last_column_entry,s1;
00707   Value m1,m2;
00708   
00709   /* Initialize all the 'Value' variables */
00710   value_init(lcm); value_init(old_lcm); value_init(gcd);
00711   value_init(last_column_entry); value_init(s1); 
00712   value_init(m1); value_init(m2);
00713 
00714   /* Compute the LCM of last column entries (denominators) of rows */
00715   value_assign(lcm,P->p[0][last_column_index]); 
00716   for(k=1;k<P->NbRows;++k) {
00717     value_assign(old_lcm,lcm);
00718     value_assign(last_column_entry,P->p[k][last_column_index]);
00719     value_gcd(gcd, lcm, last_column_entry);
00720     value_divexact(m1, last_column_entry, gcd);
00721     value_multiply(lcm,lcm,m1);
00722   }
00723   
00724   /* S[i][j] = Sum(X[i][k] * P[k][j] where Sum is extended over k = 1..nbrows*/
00725   for(i=0;i<X->NbRows;++i)
00726     for(j=0;j<P->NbColumns-1;++j) {
00727       
00728       /* Initialize s1 to zero. */
00729       value_set_si(s1,0);
00730       for(k=0;k<P->NbRows;++k) {
00731         
00732         /* If the LCM of last column entries is one, simply add the products */
00733         if(value_one_p(lcm)) {
00734           value_addmul(s1, X->p[i][k], P->p[k][j]);
00735         }  
00736         
00737         /* Numerator (num) and denominator (denom) of S[i][j] is given by :- */
00738         /* numerator  = Sum(X[i][k]*P[k][j]*lcm/P[k][last_column_index]) and */
00739         /* denominator= lcm where Sum is extended over k = 1..nbrows.        */
00740         else {
00741           value_multiply(m1,X->p[i][k],P->p[k][j]);
00742           value_division(m2,lcm,P->p[k][last_column_index]);
00743           value_addmul(s1, m1, m2);
00744         }
00745       } 
00746       value_assign(S->p[i][j],s1);
00747     }
00748   
00749   for(i=0;i<S->NbRows;++i) {
00750     value_assign(S->p[i][last_column_index],lcm);
00751 
00752     /* Normalize the rows so that last element >=0 */
00753     Vector_Normalize_Positive(&S->p[i][0],S->NbColumns,S->NbColumns-1);
00754   }
00755   
00756   /* Clear all the 'Value' variables */
00757   value_clear(lcm); value_clear(old_lcm); value_clear(gcd);
00758   value_clear(last_column_entry); value_clear(s1); 
00759   value_clear(m1); value_clear(m2);
00760  
00761   return;
00762 } /* rat_prodmat */
00763 
00764 /*
00765  * Given a matrix 'Mat' and vector 'p1', compute the matrix-vector product 
00766  * and store the result in vector 'p2'. 
00767  */
00768 void Matrix_Vector_Product(Matrix *Mat,Value *p1,Value *p2) {
00769 
00770   int NbRows, NbColumns, i, j;
00771   Value **cm, *q, *cp1, *cp2;
00772   
00773   NbRows=Mat->NbRows;
00774   NbColumns=Mat->NbColumns;
00775   
00776   cm = Mat->p;
00777   cp2 = p2;
00778   for(i=0;i<NbRows;i++) {
00779     q = *cm++;
00780     cp1 = p1;
00781     value_multiply(*cp2,*q,*cp1);
00782     q++;
00783     cp1++;
00784     
00785     /* *cp2 = *q++ * *cp1++ */
00786     for(j=1;j<NbColumns;j++) {
00787       value_addmul(*cp2, *q, *cp1);
00788       q++;
00789       cp1++;
00790     }
00791     cp2++;
00792   }
00793   return;
00794 } /* Matrix_Vector_Product */
00795 
00796 /*
00797  * Given a vector 'p1' and a matrix 'Mat', compute the vector-matrix product 
00798  * and store the result in vector 'p2'
00799  */
00800 void Vector_Matrix_Product(Value *p1,Matrix *Mat,Value *p2) {
00801   
00802   int NbRows, NbColumns, i, j;
00803   Value **cm, *cp1, *cp2;
00804   
00805   NbRows=Mat->NbRows;
00806   NbColumns=Mat->NbColumns;
00807   cp2 = p2;
00808   cm  = Mat->p;
00809   for (j=0;j<NbColumns;j++) {
00810     cp1 = p1;
00811     value_multiply(*cp2,*(*cm+j),*cp1);
00812     cp1++;
00813     
00814     /* *cp2= *(*cm+j) * *cp1++; */
00815     for (i=1;i<NbRows;i++) {
00816       value_addmul(*cp2, *(*(cm+i)+j), *cp1);
00817       cp1++;
00818     }
00819     cp2++;
00820   }
00821   return;
00822 } /* Vector_Matrix_Product */
00823 
00824 /* 
00825  * Given matrices 'Mat1' and 'Mat2', compute the matrix product and store in 
00826  * matrix 'Mat3' 
00827  */
00828 void Matrix_Product(Matrix *Mat1,Matrix *Mat2,Matrix *Mat3) {
00829   
00830   int Size, i, j, k;
00831   unsigned NbRows, NbColumns;
00832   Value **q1, **q2, *p1, *p3,sum;
00833   
00834   NbRows    = Mat1->NbRows;
00835   NbColumns = Mat2->NbColumns;
00836   
00837   Size      = Mat1->NbColumns;
00838   if(Mat2->NbRows!=Size||Mat3->NbRows!=NbRows||Mat3->NbColumns!=NbColumns) {
00839     fprintf(stderr, "? Matrix_Product : incompatable matrix dimension\n");
00840     return;
00841   }     
00842   value_init(sum); 
00843   p3 = Mat3->p_Init;
00844   q1 = Mat1->p;
00845   q2 = Mat2->p;
00846   
00847   /* Mat3[i][j] = Sum(Mat1[i][k]*Mat2[k][j] where sum is over k = 1..nbrows */
00848   for (i=0;i<NbRows;i++) {
00849     for (j=0;j<NbColumns;j++) {
00850       p1 = *(q1+i);
00851       value_set_si(sum,0);
00852       for (k=0;k<Size;k++) {
00853         value_addmul(sum, *p1, *(*(q2+k)+j));
00854         p1++;
00855       }
00856       value_assign(*p3,sum);
00857       p3++;
00858     }
00859   }
00860   value_clear(sum); 
00861   return;
00862 } /* Matrix_Product */
00863   
00864 /*
00865  * Given a rational matrix 'Mat'(k x k), compute its inverse rational matrix 
00866  * 'MatInv' k x k.
00867  * The output is 1,
00868  * if 'Mat' is non-singular (invertible), otherwise the output is 0. Note:: 
00869  * (1) Matrix 'Mat' is modified during the inverse operation.
00870  * (2) Matrix 'MatInv' must be preallocated before passing into this function.
00871  */
00872 int Matrix_Inverse(Matrix *Mat,Matrix *MatInv ) {
00873   
00874   int i, k, j, c;
00875   Value x, gcd, piv;
00876   Value m1,m2;
00877   Value *den;
00878   
00879   if(Mat->NbRows != Mat->NbColumns) {
00880    fprintf(stderr,"Trying to invert a non-square matrix !\n");
00881     return 0;
00882   }
00883   
00884   /* Initialize all the 'Value' variables */
00885   value_init(x);  value_init(gcd); value_init(piv);
00886   value_init(m1); value_init(m2);
00887 
00888   k = Mat->NbRows; 
00889 
00890   /* Initialise MatInv */
00891   Vector_Set(MatInv->p[0],0,k*k);
00892 
00893   /* Initialize 'MatInv' to Identity matrix form. Each diagonal entry is set*/
00894   /* to 1. Last column of each row (denominator of each entry in a row) is  */
00895   /* also set to 1.                                                         */ 
00896   for(i=0;i<k;++i) {
00897     value_set_si(MatInv->p[i][i],1);    
00898     /* value_set_si(MatInv->p[i][k],1); /* denum */
00899   }  
00900   /* Apply Gauss-Jordan elimination method on the two matrices 'Mat' and  */
00901   /* 'MatInv' in parallel.                                                */
00902   for(i=0;i<k;++i) {
00903     
00904     /* Check if the diagonal entry (new pivot) is non-zero or not */
00905     if(value_zero_p(Mat->p[i][i])) {    
00906       
00907       /* Search for a non-zero pivot down the column(i) */
00908       for(j=i;j<k;++j)      
00909         if(value_notzero_p(Mat->p[j][i]))
00910           break;
00911       
00912       /* If no non-zero pivot is found, the matrix 'Mat' is non-invertible */
00913       /* Return 0.                                                         */
00914       if(j==k) {
00915         
00916         /* Clear all the 'Value' variables */
00917         value_clear(x);  value_clear(gcd); value_clear(piv);
00918         value_clear(m1); value_clear(m2);
00919         return 0;
00920       } 
00921       
00922       /* Exchange the rows, row(i) and row(j) so that the diagonal element */
00923       /* Mat->p[i][i] (pivot) is non-zero. Repeat the same operations on    */
00924       /* matrix 'MatInv'.                                                   */
00925       for(c=0;c<k;++c) {
00926 
00927         /* Interchange rows, row(i) and row(j) of matrix 'Mat'    */
00928         value_assign(x,Mat->p[j][c]);
00929         value_assign(Mat->p[j][c],Mat->p[i][c]);
00930         value_assign(Mat->p[i][c],x);
00931         
00932         /* Interchange rows, row(i) and row(j) of matrix 'MatInv' */
00933         value_assign(x,MatInv->p[j][c]);
00934         value_assign(MatInv->p[j][c],MatInv->p[i][c]);
00935         value_assign(MatInv->p[i][c],x);
00936       }
00937     }
00938     
00939     /* Make all the entries in column(i) of matrix 'Mat' zero except the */
00940     /* diagonal entry. Repeat the same sequence of operations on matrix  */
00941     /* 'MatInv'.                                                         */
00942     for(j=0;j<k;++j) {
00943       if (j==i) continue;                /* Skip the pivot */
00944       value_assign(x,Mat->p[j][i]);
00945       if(value_notzero_p(x)) {
00946         value_assign(piv,Mat->p[i][i]);
00947         value_gcd(gcd, x, piv);
00948         if (value_notone_p(gcd) ) {
00949           value_divexact(x, x, gcd);
00950           value_divexact(piv, piv, gcd);
00951         }
00952         for(c=((j>i)?i:0);c<k;++c) {
00953           value_multiply(m1,piv,Mat->p[j][c]);
00954           value_multiply(m2,x,Mat->p[i][c]);
00955           value_subtract(Mat->p[j][c],m1,m2); 
00956         }
00957         for(c=0;c<k;++c) {
00958           value_multiply(m1,piv,MatInv->p[j][c]);
00959           value_multiply(m2,x,MatInv->p[i][c]);
00960           value_subtract(MatInv->p[j][c],m1,m2);
00961         }
00962               
00963         /* Simplify row(j) of the two matrices 'Mat' and 'MatInv' by */
00964         /* dividing the rows with the common GCD.                     */
00965         Vector_Gcd(&MatInv->p[j][0],k,&m1);
00966         Vector_Gcd(&Mat->p[j][0],k,&m2);
00967         value_gcd(gcd, m1, m2);
00968         if(value_notone_p(gcd)) {
00969           for(c=0;c<k;++c) {
00970             value_divexact(Mat->p[j][c], Mat->p[j][c], gcd);
00971             value_divexact(MatInv->p[j][c], MatInv->p[j][c], gcd);
00972           }
00973         }
00974       }
00975     }
00976   }
00977   
00978   /* Find common denom for each row */ 
00979    den = (Value *)malloc(k*sizeof(Value));
00980    value_set_si(x,1);
00981    for(j=0 ; j<k ; ++j) {
00982      value_init(den[j]);
00983      value_assign(den[j],Mat->p[j][j]);
00984      
00985      /* gcd is always positive */
00986      Vector_Gcd(&MatInv->p[j][0],k,&gcd);
00987      value_gcd(gcd, gcd, den[j]);
00988      if (value_neg_p(den[j])) 
00989        value_oppose(gcd,gcd); /* make denominator positive */
00990      if (value_notone_p(gcd)) {
00991        for (c=0; c<k; c++) 
00992          value_divexact(MatInv->p[j][c], MatInv->p[j][c], gcd); /* normalize */
00993        value_divexact(den[j], den[j], gcd);
00994      }  
00995      value_gcd(gcd, x, den[j]);
00996      value_divexact(m1, den[j], gcd);
00997      value_multiply(x,x,m1);
00998    }
00999    if (value_notone_p(x)) 
01000      for(j=0 ; j<k ; ++j) {       
01001        for (c=0; c<k; c++) {
01002          value_division(m1,x,den[j]);
01003          value_multiply(MatInv->p[j][c],MatInv->p[j][c],m1);  /* normalize */
01004        }
01005      }
01006 
01007    /* Clear all the 'Value' variables */
01008    for(j=0 ; j<k ; ++j) {
01009      value_clear(den[j]);
01010    }  
01011    value_clear(x);  value_clear(gcd); value_clear(piv);
01012    value_clear(m1); value_clear(m2);
01013    free(den);
01014    
01015    return 1;
01016 } /* Matrix_Inverse */
01017 
01018 
01019 
01020 
01021 
01022 
01023 
01024 

Generated on Thu Sep 4 15:28:58 2008 for polylib by doxygen 1.3.5