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