polylib 5.22.8
matrix.c
Go to the documentation of this file.
1/* matrix.c
2 COPYRIGHT
3 Both this software and its documentation are
4
5 Copyright 1993 by IRISA /Universite de Rennes I -
6 France, Copyright 1995,1996 by BYU, Provo, Utah
7 all rights reserved.
8
9*/
10
11#include <ctype.h>
12#include <polylib/polylib.h>
13#include <stdio.h>
14#include <stdlib.h>
15#include <string.h>
16
17#ifdef mac_os
18#define abs __abs
19#endif
20
21/*
22 * Allocate space for matrix dimensioned by 'NbRows X NbColumns'.
23 */
24Matrix *Matrix_Alloc(unsigned NbRows, unsigned NbColumns) {
25
26 Matrix *Mat;
27 Value *p, **q;
28 int i;
29
30 Mat = (Matrix *)malloc(sizeof(Matrix));
31 if (!Mat) {
32 errormsg1("Matrix_Alloc", "outofmem", "out of memory space");
33 return 0;
34 }
35 Mat->NbRows = NbRows;
36 Mat->NbColumns = NbColumns;
37 if (NbRows == 0 || NbColumns == 0) {
38 Mat->p = (Value **)0;
39 Mat->p_Init = (Value *)0;
40 Mat->p_Init_size = 0;
41 } else {
42 q = (Value **)malloc(NbRows * sizeof(*q));
43 if (!q) {
44 free(Mat);
45 errormsg1("Matrix_Alloc", "outofmem", "out of memory space");
46 return 0;
47 }
48 p = value_alloc(NbRows * NbColumns, &Mat->p_Init_size);
49 if (!p) {
50 free(q);
51 free(Mat);
52 errormsg1("Matrix_Alloc", "outofmem", "out of memory space");
53 return 0;
54 }
55 Mat->p = q;
56 Mat->p_Init = p;
57 for (i = 0; i < NbRows; i++) {
58 *q++ = p;
59 p += NbColumns;
60 }
61 }
62 p = NULL;
63 q = NULL;
64
65 return Mat;
66} /* Matrix_Alloc */
67
68/*
69 * Free the memory space occupied by Matrix 'Mat'
70 */
71void Matrix_Free(Matrix *Mat) {
72 if (Mat->p_Init)
73 value_free(Mat->p_Init, Mat->p_Init_size);
74
75 if (Mat->p)
76 free(Mat->p);
77 free(Mat);
78
79} /* Matrix_Free */
80
81void Matrix_Extend(Matrix *Mat, unsigned NbRows) {
82 Value *p, **q;
83 int i;
84
85 q = (Value **)realloc(Mat->p, NbRows * sizeof(*q));
86 if (!q) {
87 errormsg1("Matrix_Extend", "outofmem", "out of memory space");
88 return;
89 }
90 Mat->p = q;
91 if (Mat->p_Init_size < NbRows * Mat->NbColumns) {
92 p = (Value *)realloc(Mat->p_Init, NbRows * Mat->NbColumns * sizeof(Value));
93 if (!p) {
94 errormsg1("Matrix_Extend", "outofmem", "out of memory space");
95 return;
96 }
97 Mat->p_Init = p;
98 Vector_Set(Mat->p_Init + Mat->NbRows * Mat->NbColumns, 0,
99 Mat->p_Init_size - Mat->NbRows * Mat->NbColumns);
100 for (i = Mat->p_Init_size; i < Mat->NbColumns * NbRows; ++i)
101 value_init(Mat->p_Init[i]);
102 Mat->p_Init_size = Mat->NbColumns * NbRows;
103 } else
104 Vector_Set(Mat->p_Init + Mat->NbRows * Mat->NbColumns, 0,
105 (NbRows - Mat->NbRows) * Mat->NbColumns);
106 for (i = 0; i < NbRows; i++) {
107 Mat->p[i] = Mat->p_Init + (i * Mat->NbColumns);
108 }
109 Mat->NbRows = NbRows;
110}
111
112/*
113 * Print the contents of the Matrix 'Mat'
114 */
115void Matrix_Print(FILE *Dst, const char *Format, Matrix *Mat) {
116 Value *p;
117 int i, j;
118 unsigned NbRows, NbColumns;
119
120 fprintf(Dst, "%d %d\n", NbRows = Mat->NbRows, NbColumns = Mat->NbColumns);
121 if (NbColumns == 0) {
122 fprintf(Dst, "\n");
123 return;
124 }
125 for (i = 0; i < NbRows; i++) {
126 p = *(Mat->p + i);
127 for (j = 0; j < NbColumns; j++) {
128 if (!Format) {
129 value_print(Dst, " " P_VALUE_FMT " ", *p++);
130 } else {
131 value_print(Dst, Format, *p++);
132 }
133 }
134 fprintf(Dst, "\n");
135 }
136} /* Matrix_Print */
137
138/*
139 * Read the contents of the Matrix 'Mat' from standard input
140 */
142 return Matrix_Read_InputFile(Mat, stdin);
143} /* Matrix_Read_Input */
144
145/*
146 * Read the contents of the Matrix 'Mat' from file open in fp
147 */
149
150 Value *p;
151 int i, j, n;
152 char *c, s[1024];
153
154 p = Mat->p_Init;
155 for (i = 0; i < Mat->NbRows; i++) {
156 do {
157 c = fgets(s, 1024, fp);
158 if (!c)
159 break;
160 /* jump to the first non space char */
161 while (isspace(*c) && *c != '\n' && *c)
162 ++c;
163 } while (*c == '#' ||
164 *c == '\n'); /* continue if the line is empty or is a comment */
165 if (!c) {
166 errormsg1("Matrix_Read", "baddim", "not enough rows");
167 return (NULL);
168 }
169
170 if (c - s >= 1023) {
171 /* skip to EOL and ignore the rest of the line if longer than 1024 */
172 errormsg1(
173 "Matrix_Read", "warning",
174 "line longer than 1024 chars (ignored remaining chars till EOL)");
175 do
176 n = getchar();
177 while (n != EOF && n != '\n');
178 }
179
180 for (j = 0; j < Mat->NbColumns; j++) {
181 char *z;
182 if (*c == '\n' || *c == '#' || *c == '\0') {
183 errormsg1("Matrix_Read", "baddim", "not enough columns");
184 return (NULL);
185 }
186 /* go the the next space or end */
187 for (z = c; *z; z++) {
188 if (*z == '\n' || *z == '#' || isspace(*z))
189 break;
190 }
191 if (*z)
192 *z = '\0';
193 else
194 z--; /* hit EOL :/ go back one char */
195 value_read(*(p++), c);
196 /* go point to the next non space char */
197 c = z + 1;
198 while (isspace(*c))
199 c++;
200 }
201 }
202 return (Mat);
203} /* Matrix_Read_InputFile */
204
205/*
206 * Read the contents of the matrix 'Mat' from standard input.
207 * a '#' is a comment (till EOL)
208 */
209Matrix *Matrix_Read(void) { return Matrix_ReadFile(stdin); } /* Matrix_Read */
210
211/*
212 * Read the contents of the matrix 'Mat' from file open in fp.
213 * a '#' is a comment (till EOL)
214 */
216 Matrix *Mat;
217 unsigned NbRows, NbColumns;
218 char s[1024];
219
220 if (fgets(s, 1024, fp) == NULL)
221 return NULL;
222 while ((*s == '#' || *s == '\n') ||
223 (sscanf(s, "%d %d", &NbRows, &NbColumns) < 2)) {
224 if (fgets(s, 1024, fp) == NULL)
225 return NULL;
226 }
227 Mat = Matrix_Alloc(NbRows, NbColumns);
228 if (!Mat) {
229 errormsg1("Matrix_Read", "outofmem", "out of memory space");
230 return (NULL);
231 }
232 if (!Matrix_Read_Input(Mat)) {
233 Matrix_Free(Mat);
234 Mat = NULL;
235 }
236 return Mat;
237} /* Matrix_ReadFile */
238
239/*
240 * Basic hermite engine
241 */
242static int hermite(Matrix *H, Matrix *U, Matrix *Q) {
243
244 int nc, nr, i, j, k, rank, reduced, pivotrow;
245 Value pivot, x, aux;
246 Value *temp1, *temp2;
247
248 /* T -1 T */
249 /* Computes form: A = Q H and U A = H and U = Q */
250
251 if (!H) {
252 errormsg1("Domlib", "nullH", "hermite: ? Null H");
253 return -1;
254 }
255 nc = H->NbColumns;
256 nr = H->NbRows;
257 temp1 = (Value *)malloc(nc * sizeof(Value));
258 temp2 = (Value *)malloc(nr * sizeof(Value));
259 if (!temp1 || !temp2) {
260 errormsg1("Domlib", "outofmem", "out of memory space");
261 return -1;
262 }
263
264 /* Initialize all the 'Value' variables */
265 value_init(pivot);
266 value_init(x);
267 value_init(aux);
268 for (i = 0; i < nc; i++)
269 value_init(temp1[i]);
270 for (i = 0; i < nr; i++)
271 value_init(temp2[i]);
272
273#ifdef DEBUG
274 fprintf(stderr, "Start -----------\n");
275 Matrix_Print(stderr, 0, H);
276#endif
277 for (k = 0, rank = 0; k < nc && rank < nr; k = k + 1) {
278 reduced = 1; /* go through loop the first time */
279#ifdef DEBUG
280 fprintf(stderr, "Working on col %d. Rank=%d ----------\n", k + 1,
281 rank + 1);
282#endif
283 while (reduced) {
284 reduced = 0;
285
286 /* 1. find pivot row */
287 value_absolute(pivot, H->p[rank][k]);
288
289 /* the kth-diagonal element */
290 pivotrow = rank;
291
292 /* find the row i>rank with smallest nonzero element in col k */
293 for (i = rank + 1; i < nr; i++) {
294 value_absolute(x, H->p[i][k]);
295 if (value_notzero_p(x) && (value_lt(x, pivot) || value_zero_p(pivot))) {
296 value_assign(pivot, x);
297 pivotrow = i;
298 }
299 }
300
301 /* 2. Bring pivot to diagonal (exchange rows pivotrow and rank) */
302 if (pivotrow != rank) {
303 Vector_Exchange(H->p[pivotrow], H->p[rank], nc);
304 if (U)
305 Vector_Exchange(U->p[pivotrow], U->p[rank], nr);
306 if (Q)
307 Vector_Exchange(Q->p[pivotrow], Q->p[rank], nr);
308
309#ifdef DEBUG
310 fprintf(stderr, "Exchange rows %d and %d -----------\n", rank + 1,
311 pivotrow + 1);
312 Matrix_Print(stderr, 0, H);
313#endif
314 }
315 value_assign(pivot, H->p[rank][k]); /* actual ( no abs() ) pivot */
316
317 /* 3. Invert the row 'rank' if pivot is negative */
318 if (value_neg_p(pivot)) {
319 value_oppose(pivot, pivot); /* pivot = -pivot */
320 for (j = 0; j < nc; j++)
321 value_oppose(H->p[rank][j], H->p[rank][j]);
322
323 /* H->p[rank][j] = -(H->p[rank][j]); */
324 if (U)
325 for (j = 0; j < nr; j++)
326 value_oppose(U->p[rank][j], U->p[rank][j]);
327
328 /* U->p[rank][j] = -(U->p[rank][j]); */
329 if (Q)
330 for (j = 0; j < nr; j++)
331 value_oppose(Q->p[rank][j], Q->p[rank][j]);
332
333 /* Q->p[rank][j] = -(Q->p[rank][j]); */
334#ifdef DEBUG
335 fprintf(stderr, "Negate row %d -----------\n", rank + 1);
336 Matrix_Print(stderr, 0, H);
337#endif
338 }
339 if (value_notzero_p(pivot)) {
340
341 /* 4. Reduce the column modulo the pivot */
342 /* This eventually zeros out everything below the */
343 /* diagonal and produces an upper triangular matrix */
344
345 for (i = rank + 1; i < nr; i++) {
346 value_assign(x, H->p[i][k]);
347 if (value_notzero_p(x)) {
348 value_modulus(aux, x, pivot);
349
350 /* floor[integer division] (corrected for neg x) */
351 if (value_neg_p(x) && value_notzero_p(aux)) {
352
353 /* x=(x/pivot)-1; */
354 value_division(x, x, pivot);
355 value_decrement(x, x);
356 } else
357 value_division(x, x, pivot);
358 for (j = 0; j < nc; j++) {
359 value_multiply(aux, x, H->p[rank][j]);
360 value_subtract(H->p[i][j], H->p[i][j], aux);
361 }
362
363 /* U->p[i][j] -= (x * U->p[rank][j]); */
364 if (U)
365 for (j = 0; j < nr; j++) {
366 value_multiply(aux, x, U->p[rank][j]);
367 value_subtract(U->p[i][j], U->p[i][j], aux);
368 }
369
370 /* Q->p[rank][j] += (x * Q->p[i][j]); */
371 if (Q)
372 for (j = 0; j < nr; j++) {
373 value_addmul(Q->p[rank][j], x, Q->p[i][j]);
374 }
375 reduced = 1;
376
377#ifdef DEBUG
378 fprintf(stderr, "row %d = row %d - %d row %d -----------\n", i + 1,
379 i + 1, x, rank + 1);
380 Matrix_Print(stderr, 0, H);
381#endif
382
383 } /* if (x) */
384 } /* for (i) */
385 } /* if (pivot != 0) */
386 } /* while (reduced) */
387
388 /* Last finish up this column */
389 /* 5. Make pivot column positive (above pivot row) */
390 /* x should be zero for i>k */
391
392 if (value_notzero_p(pivot)) {
393 for (i = 0; i < rank; i++) {
394 value_assign(x, H->p[i][k]);
395 if (value_notzero_p(x)) {
396 value_modulus(aux, x, pivot);
397
398 /* floor[integer division] (corrected for neg x) */
399 if (value_neg_p(x) && value_notzero_p(aux)) {
400 value_division(x, x, pivot);
401 value_decrement(x, x);
402
403 /* x=(x/pivot)-1; */
404 } else
405 value_division(x, x, pivot);
406
407 /* H->p[i][j] -= x * H->p[rank][j]; */
408 for (j = 0; j < nc; j++) {
409 value_multiply(aux, x, H->p[rank][j]);
410 value_subtract(H->p[i][j], H->p[i][j], aux);
411 }
412
413 /* U->p[i][j] -= x * U->p[rank][j]; */
414 if (U)
415 for (j = 0; j < nr; j++) {
416 value_multiply(aux, x, U->p[rank][j]);
417 value_subtract(U->p[i][j], U->p[i][j], aux);
418 }
419
420 /* Q->p[rank][j] += x * Q->p[i][j]; */
421 if (Q)
422 for (j = 0; j < nr; j++) {
423 value_addmul(Q->p[rank][j], x, Q->p[i][j]);
424 }
425#ifdef DEBUG
426 fprintf(stderr, "row %d = row %d - %d row %d -----------\n", i + 1,
427 i + 1, x, rank + 1);
428 Matrix_Print(stderr, 0, H);
429#endif
430 } /* if (x) */
431 } /* for (i) */
432 rank++;
433 } /* if (pivot!=0) */
434 } /* for (k) */
435
436 /* Clear all the 'Value' variables */
437 value_clear(pivot);
438 value_clear(x);
439 value_clear(aux);
440 for (i = 0; i < nc; i++)
441 value_clear(temp1[i]);
442 for (i = 0; i < nr; i++)
443 value_clear(temp2[i]);
444 free(temp2);
445 free(temp1);
446 return rank;
447} /* Hermite */
448
449void right_hermite(Matrix *A, Matrix **Hp, Matrix **Up, Matrix **Qp) {
450
451 Matrix *H, *Q, *U;
452 int i, j, nr, nc;
453 Value tmp;
454
455 /* Computes form: A = QH , UA = H */
456 nc = A->NbColumns;
457 nr = A->NbRows;
458
459 /* H = A */
460 *Hp = H = Matrix_Alloc(nr, nc);
461 if (!H) {
462 errormsg1("DomRightHermite", "outofmem", "out of memory space");
463 return;
464 }
465
466 /* Initialize all the 'Value' variables */
467 value_init(tmp);
468
469 Vector_Copy(A->p_Init, H->p_Init, nr * nc);
470
471 /* U = I */
472 if (Up) {
473 *Up = U = Matrix_Alloc(nr, nr);
474 if (!U) {
475 errormsg1("DomRightHermite", "outofmem", "out of memory space");
476 value_clear(tmp);
477 return;
478 }
479 Vector_Set(U->p_Init, 0, nr * nr); /* zero's */
480 for (i = 0; i < nr; i++) /* with diagonal of 1's */
481 value_set_si(U->p[i][i], 1);
482 } else
483 U = (Matrix *)0;
484
485 /* Q = I */
486 /* Actually I compute Q transpose... its easier */
487 if (Qp) {
488 *Qp = Q = Matrix_Alloc(nr, nr);
489 if (!Q) {
490 errormsg1("DomRightHermite", "outofmem", "out of memory space");
491 value_clear(tmp);
492 return;
493 }
494 Vector_Set(Q->p_Init, 0, nr * nr); /* zero's */
495 for (i = 0; i < nr; i++) /* with diagonal of 1's */
496 value_set_si(Q->p[i][i], 1);
497 } else
498 Q = (Matrix *)0;
499
500 hermite(H, U, Q);
501
502 /* Q is returned transposed */
503 /* Transpose Q */
504 if (Q) {
505 for (i = 0; i < nr; i++) {
506 for (j = i + 1; j < nr; j++) {
507 value_assign(tmp, Q->p[i][j]);
508 value_assign(Q->p[i][j], Q->p[j][i]);
509 value_assign(Q->p[j][i], tmp);
510 }
511 }
512 }
513 value_clear(tmp);
514 return;
515} /* right_hermite */
516
517void left_hermite(Matrix *A, Matrix **Hp, Matrix **Qp, Matrix **Up) {
518
519 Matrix *H, *HT, *Q, *U;
520 int i, j, nc, nr;
521 Value tmp;
522
523 /* Computes left form: A = HQ , AU = H ,
524 T T T T T T
525 using right form A = Q H , U A = H */
526
527 nr = A->NbRows;
528 nc = A->NbColumns;
529
530 /* HT = A transpose */
531 HT = Matrix_Alloc(nc, nr);
532 if (!HT) {
533 errormsg1("DomLeftHermite", "outofmem", "out of memory space");
534 return;
535 }
536 value_init(tmp);
537 for (i = 0; i < nr; i++)
538 for (j = 0; j < nc; j++)
539 value_assign(HT->p[j][i], A->p[i][j]);
540
541 /* U = I */
542 if (Up) {
543 *Up = U = Matrix_Alloc(nc, nc);
544 if (!U) {
545 errormsg1("DomLeftHermite", "outofmem", "out of memory space");
546 value_clear(tmp);
547 return;
548 }
549 Vector_Set(U->p_Init, 0, nc * nc); /* zero's */
550 for (i = 0; i < nc; i++) /* with diagonal of 1's */
551 value_set_si(U->p[i][i], 1);
552 } else
553 U = (Matrix *)0;
554
555 /* Q = I */
556 if (Qp) {
557 *Qp = Q = Matrix_Alloc(nc, nc);
558 if (!Q) {
559 errormsg1("DomLeftHermite", "outofmem", "out of memory space");
560 value_clear(tmp);
561 return;
562 }
563 Vector_Set(Q->p_Init, 0, nc * nc); /* zero's */
564 for (i = 0; i < nc; i++) /* with diagonal of 1's */
565 value_set_si(Q->p[i][i], 1);
566 } else
567 Q = (Matrix *)0;
568 hermite(HT, U, Q);
569
570 /* H = HT transpose */
571 *Hp = H = Matrix_Alloc(nr, nc);
572 if (!H) {
573 errormsg1("DomLeftHermite", "outofmem", "out of memory space");
574 value_clear(tmp);
575 return;
576 }
577 for (i = 0; i < nr; i++)
578 for (j = 0; j < nc; j++)
579 value_assign(H->p[i][j], HT->p[j][i]);
580 Matrix_Free(HT);
581
582 /* Transpose U */
583 if (U) {
584 for (i = 0; i < nc; i++) {
585 for (j = i + 1; j < nc; j++) {
586 value_assign(tmp, U->p[i][j]);
587 value_assign(U->p[i][j], U->p[j][i]);
588 value_assign(U->p[j][i], tmp);
589 }
590 }
591 }
592 value_clear(tmp);
593} /* left_hermite */
594
595/*
596 * Given a integer matrix 'Mat'(k x k), compute its inverse rational matrix
597 * 'MatInv' k x (k+1). The last column of each row in matrix MatInv is used
598 * to store the common denominator of the entries in a row. The output is 1,
599 * if 'Mat' is non-singular (invertible), otherwise the output is 0. Note::
600 * (1) Matrix 'Mat' is modified during the inverse operation.
601 * (2) Matrix 'MatInv' must be preallocated before passing into this function.
602 */
603int MatInverse(Matrix *Mat, Matrix *MatInv) {
604
605 int i, k, j, c;
606 Value x, gcd, piv;
607 Value m1, m2;
608
609 if (Mat->NbRows != Mat->NbColumns) {
610 fprintf(stderr, "Trying to invert a non-square matrix !\n");
611 return 0;
612 }
613
614 /* Initialize all the 'Value' variables */
615 value_init(x);
616 value_init(gcd);
617 value_init(piv);
618 value_init(m1);
619 value_init(m2);
620
621 k = Mat->NbRows;
622
623 /* Initialise MatInv */
624 Vector_Set(MatInv->p[0], 0, k * (k + 1));
625
626 /* Initialize 'MatInv' to Identity matrix form. Each diagonal entry is set*/
627 /* to 1. Last column of each row (denominator of each entry in a row) is */
628 /* also set to 1. */
629 for (i = 0; i < k; ++i) {
630 value_set_si(MatInv->p[i][i], 1);
631 value_set_si(MatInv->p[i][k], 1); /* denum */
632 }
633 /* Apply Gauss-Jordan elimination method on the two matrices 'Mat' and */
634 /* 'MatInv' in parallel. */
635 for (i = 0; i < k; ++i) {
636
637 /* Check if the diagonal entry (new pivot) is non-zero or not */
638 if (value_zero_p(Mat->p[i][i])) {
639
640 /* Search for a non-zero pivot down the column(i) */
641 for (j = i; j < k; ++j)
642 if (value_notzero_p(Mat->p[j][i]))
643 break;
644
645 /* If no non-zero pivot is found, the matrix 'Mat' is non-invertible */
646 /* Return 0. */
647 if (j == k) {
648
649 /* Clear all the 'Value' variables */
650 value_clear(x);
651 value_clear(gcd);
652 value_clear(piv);
653 value_clear(m1);
654 value_clear(m2);
655 return 0;
656 }
657
658 /* Exchange the rows, row(i) and row(j) so that the diagonal element */
659 /* Mat->p[i][i] (pivot) is non-zero. Repeat the same operations on */
660 /* matrix 'MatInv'. */
661 for (c = 0; c < k; ++c) {
662
663 /* Interchange rows, row(i) and row(j) of matrix 'Mat' */
664 value_assign(x, Mat->p[j][c]);
665 value_assign(Mat->p[j][c], Mat->p[i][c]);
666 value_assign(Mat->p[i][c], x);
667
668 /* Interchange rows, row(i) and row(j) of matrix 'MatInv' */
669 value_assign(x, MatInv->p[j][c]);
670 value_assign(MatInv->p[j][c], MatInv->p[i][c]);
671 value_assign(MatInv->p[i][c], x);
672 }
673 }
674
675 /* Make all the entries in column(i) of matrix 'Mat' zero except the */
676 /* diagonal entry. Repeat the same sequence of operations on matrix */
677 /* 'MatInv'. */
678 for (j = 0; j < k; ++j) {
679 if (j == i)
680 continue; /* Skip the pivot */
681 value_assign(x, Mat->p[j][i]);
682 if (value_notzero_p(x)) {
683 value_assign(piv, Mat->p[i][i]);
684 value_gcd(gcd, x, piv);
685 if (value_notone_p(gcd)) {
686 value_divexact(x, x, gcd);
687 value_divexact(piv, piv, gcd);
688 }
689 for (c = ((j > i) ? i : 0); c < k; ++c) {
690 value_multiply(m1, piv, Mat->p[j][c]);
691 value_multiply(m2, x, Mat->p[i][c]);
692 value_subtract(Mat->p[j][c], m1, m2);
693 }
694 for (c = 0; c < k; ++c) {
695 value_multiply(m1, piv, MatInv->p[j][c]);
696 value_multiply(m2, x, MatInv->p[i][c]);
697 value_subtract(MatInv->p[j][c], m1, m2);
698 }
699
700 /* Simplify row(j) of the two matrices 'Mat' and 'MatInv' by */
701 /* dividing the rows with the common GCD. */
702 Vector_Gcd(&MatInv->p[j][0], k, &m1);
703 Vector_Gcd(&Mat->p[j][0], k, &m2);
704 value_gcd(gcd, m1, m2);
705 if (value_notone_p(gcd)) {
706 for (c = 0; c < k; ++c) {
707 value_divexact(Mat->p[j][c], Mat->p[j][c], gcd);
708 value_divexact(MatInv->p[j][c], MatInv->p[j][c], gcd);
709 }
710 }
711 }
712 }
713 }
714
715 /* Simplify every row so that 'Mat' reduces to Identity matrix. Perform */
716 /* the same sequence of operations on the matrix 'MatInv'. */
717 for (j = 0; j < k; ++j) {
718 value_assign(MatInv->p[j][k], Mat->p[j][j]);
719
720 /* Make the last column (denominator of each entry) of every row greater */
721 /* than zero. */
722 Vector_Normalize_Positive(&MatInv->p[j][0], k + 1, k);
723 }
724
725 /* Clear all the 'Value' variables */
726 value_clear(x);
727 value_clear(gcd);
728 value_clear(piv);
729 value_clear(m1);
730 value_clear(m2);
731
732 return 1;
733} /* Mat_Inverse */
734
735/*
736 * Given (m x n) integer matrix 'X' and n x (k+1) rational matrix 'P', compute
737 * the rational m x (k+1) rational matrix 'S'. The last column in each row of
738 * the rational matrices is used to store the common denominator of elements
739 * in a row.
740 */
741void rat_prodmat(Matrix *S, Matrix *X, Matrix *P) {
742
743 int i, j, k;
744 int last_column_index = P->NbColumns - 1;
745 Value lcm, gcd, last_column_entry, s1;
746 Value m1, m2;
747
748 /* Initialize all the 'Value' variables */
749 value_init(lcm);
750 value_init(gcd);
751 value_init(last_column_entry);
752 value_init(s1);
753 value_init(m1);
754 value_init(m2);
755
756 /* Compute the LCM of last column entries (denominators) of rows */
757 value_assign(lcm, P->p[0][last_column_index]);
758 for (k = 1; k < P->NbRows; ++k) {
759 value_assign(last_column_entry, P->p[k][last_column_index]);
760 value_gcd(gcd, lcm, last_column_entry);
761 value_divexact(m1, last_column_entry, gcd);
762 value_multiply(lcm, lcm, m1);
763 }
764
765 /* S[i][j] = Sum(X[i][k] * P[k][j] where Sum is extended over k = 1..nbrows*/
766 for (i = 0; i < X->NbRows; ++i)
767 for (j = 0; j < P->NbColumns - 1; ++j) {
768
769 /* Initialize s1 to zero. */
770 value_set_si(s1, 0);
771 for (k = 0; k < P->NbRows; ++k) {
772
773 /* If the LCM of last column entries is one, simply add the products */
774 if (value_one_p(lcm)) {
775 value_addmul(s1, X->p[i][k], P->p[k][j]);
776 }
777
778 /* Numerator (num) and denominator (denom) of S[i][j] is given by :- */
779 /* numerator = Sum(X[i][k]*P[k][j]*lcm/P[k][last_column_index]) and */
780 /* denominator= lcm where Sum is extended over k = 1..nbrows. */
781 else {
782 value_multiply(m1, X->p[i][k], P->p[k][j]);
783 value_division(m2, lcm, P->p[k][last_column_index]);
784 value_addmul(s1, m1, m2);
785 }
786 }
787 value_assign(S->p[i][j], s1);
788 }
789
790 for (i = 0; i < S->NbRows; ++i) {
791 value_assign(S->p[i][last_column_index], lcm);
792
793 /* Normalize the rows so that last element >=0 */
794 Vector_Normalize_Positive(&S->p[i][0], S->NbColumns, S->NbColumns - 1);
795 }
796
797 /* Clear all the 'Value' variables */
798 value_clear(lcm);
799 value_clear(gcd);
800 value_clear(last_column_entry);
801 value_clear(s1);
802 value_clear(m1);
803 value_clear(m2);
804
805 return;
806} /* rat_prodmat */
807
808/*
809 * Given a matrix 'Mat' and vector 'p1', compute the matrix-vector product
810 * and store the result in vector 'p2'.
811 */
812void Matrix_Vector_Product(Matrix *Mat, Value *p1, Value *p2) {
813
814 int NbRows, NbColumns, i, j;
815 Value **cm, *q, *cp1, *cp2;
816
817 NbRows = Mat->NbRows;
818 NbColumns = Mat->NbColumns;
819
820 cm = Mat->p;
821 cp2 = p2;
822 for (i = 0; i < NbRows; i++) {
823 q = *cm++;
824 cp1 = p1;
825 value_multiply(*cp2, *q, *cp1);
826 q++;
827 cp1++;
828
829 /* *cp2 = *q++ * *cp1++ */
830 for (j = 1; j < NbColumns; j++) {
831 value_addmul(*cp2, *q, *cp1);
832 q++;
833 cp1++;
834 }
835 cp2++;
836 }
837 return;
838} /* Matrix_Vector_Product */
839
840/*
841 * Given a vector 'p1' and a matrix 'Mat', compute the vector-matrix product
842 * and store the result in vector 'p2'
843 */
844void Vector_Matrix_Product(Value *p1, Matrix *Mat, Value *p2) {
845
846 int NbRows, NbColumns, i, j;
847 Value **cm, *cp1, *cp2;
848
849 NbRows = Mat->NbRows;
850 NbColumns = Mat->NbColumns;
851 cp2 = p2;
852 cm = Mat->p;
853 for (j = 0; j < NbColumns; j++) {
854 cp1 = p1;
855 value_multiply(*cp2, *(*cm + j), *cp1);
856 cp1++;
857
858 /* *cp2= *(*cm+j) * *cp1++; */
859 for (i = 1; i < NbRows; i++) {
860 value_addmul(*cp2, *(*(cm + i) + j), *cp1);
861 cp1++;
862 }
863 cp2++;
864 }
865 return;
866} /* Vector_Matrix_Product */
867
868/*
869 * Given matrices 'Mat1' and 'Mat2', compute the matrix product and store in
870 * matrix 'Mat3'
871 */
872void Matrix_Product(Matrix *Mat1, Matrix *Mat2, Matrix *Mat3) {
873
874 int Size, i, j, k;
875 unsigned NbRows, NbColumns;
876 Value **q1, **q2, *p1, *p3, sum;
877
878 NbRows = Mat1->NbRows;
879 NbColumns = Mat2->NbColumns;
880
881 Size = Mat1->NbColumns;
882 if (Mat2->NbRows != Size || Mat3->NbRows != NbRows ||
883 Mat3->NbColumns != NbColumns) {
884 fprintf(stderr, "? Matrix_Product : incompatible matrix dimension\n");
885 return;
886 }
887 value_init(sum);
888 p3 = Mat3->p_Init;
889 q1 = Mat1->p;
890 q2 = Mat2->p;
891
892 /* Mat3[i][j] = Sum(Mat1[i][k]*Mat2[k][j] where sum is over k = 1..nbrows */
893 for (i = 0; i < NbRows; i++) {
894 for (j = 0; j < NbColumns; j++) {
895 p1 = *(q1 + i);
896 value_set_si(sum, 0);
897 for (k = 0; k < Size; k++) {
898 value_addmul(sum, *p1, *(*(q2 + k) + j));
899 p1++;
900 }
901 value_assign(*p3, sum);
902 p3++;
903 }
904 }
905 value_clear(sum);
906 return;
907} /* Matrix_Product */
908
909/*
910 * Given a rational matrix 'Mat'(k x k), compute its inverse rational matrix
911 * 'MatInv' k x k.
912 * The output is 1,
913 * if 'Mat' is non-singular (invertible), otherwise the output is 0. Note::
914 * (1) Matrix 'Mat' is modified during the inverse operation.
915 * (2) Matrix 'MatInv' must be preallocated before passing into this function.
916 */
917int Matrix_Inverse(Matrix *Mat, Matrix *MatInv) {
918
919 int i, k, j, c;
920 Value x, gcd, piv;
921 Value m1, m2;
922 Value *den;
923
924 if (Mat->NbRows != Mat->NbColumns) {
925 fprintf(stderr, "Trying to invert a non-square matrix !\n");
926 return 0;
927 }
928
929 /* Initialize all the 'Value' variables */
930 value_init(x);
931 value_init(gcd);
932 value_init(piv);
933 value_init(m1);
934 value_init(m2);
935
936 k = Mat->NbRows;
937
938 /* Initialise MatInv */
939 Vector_Set(MatInv->p[0], 0, k * k);
940
941 /* Initialize 'MatInv' to Identity matrix form. Each diagonal entry is set*/
942 /* to 1. Last column of each row (denominator of each entry in a row) is */
943 /* also set to 1. */
944 for (i = 0; i < k; ++i) {
945 value_set_si(MatInv->p[i][i], 1);
946 /* value_set_si(MatInv->p[i][k],1); // denum */
947 }
948 /* Apply Gauss-Jordan elimination method on the two matrices 'Mat' and */
949 /* 'MatInv' in parallel. */
950 for (i = 0; i < k; ++i) {
951
952 /* Check if the diagonal entry (new pivot) is non-zero or not */
953 if (value_zero_p(Mat->p[i][i])) {
954
955 /* Search for a non-zero pivot down the column(i) */
956 for (j = i; j < k; ++j)
957 if (value_notzero_p(Mat->p[j][i]))
958 break;
959
960 /* If no non-zero pivot is found, the matrix 'Mat' is non-invertible */
961 /* Return 0. */
962 if (j == k) {
963
964 /* Clear all the 'Value' variables */
965 value_clear(x);
966 value_clear(gcd);
967 value_clear(piv);
968 value_clear(m1);
969 value_clear(m2);
970 return 0;
971 }
972
973 /* Exchange the rows, row(i) and row(j) so that the diagonal element */
974 /* Mat->p[i][i] (pivot) is non-zero. Repeat the same operations on */
975 /* matrix 'MatInv'. */
976 for (c = 0; c < k; ++c) {
977
978 /* Interchange rows, row(i) and row(j) of matrix 'Mat' */
979 value_assign(x, Mat->p[j][c]);
980 value_assign(Mat->p[j][c], Mat->p[i][c]);
981 value_assign(Mat->p[i][c], x);
982
983 /* Interchange rows, row(i) and row(j) of matrix 'MatInv' */
984 value_assign(x, MatInv->p[j][c]);
985 value_assign(MatInv->p[j][c], MatInv->p[i][c]);
986 value_assign(MatInv->p[i][c], x);
987 }
988 }
989
990 /* Make all the entries in column(i) of matrix 'Mat' zero except the */
991 /* diagonal entry. Repeat the same sequence of operations on matrix */
992 /* 'MatInv'. */
993 for (j = 0; j < k; ++j) {
994 if (j == i)
995 continue; /* Skip the pivot */
996 value_assign(x, Mat->p[j][i]);
997 if (value_notzero_p(x)) {
998 value_assign(piv, Mat->p[i][i]);
999 value_gcd(gcd, x, piv);
1000 if (value_notone_p(gcd)) {
1001 value_divexact(x, x, gcd);
1002 value_divexact(piv, piv, gcd);
1003 }
1004 for (c = ((j > i) ? i : 0); c < k; ++c) {
1005 value_multiply(m1, piv, Mat->p[j][c]);
1006 value_multiply(m2, x, Mat->p[i][c]);
1007 value_subtract(Mat->p[j][c], m1, m2);
1008 }
1009 for (c = 0; c < k; ++c) {
1010 value_multiply(m1, piv, MatInv->p[j][c]);
1011 value_multiply(m2, x, MatInv->p[i][c]);
1012 value_subtract(MatInv->p[j][c], m1, m2);
1013 }
1014
1015 /* Simplify row(j) of the two matrices 'Mat' and 'MatInv' by */
1016 /* dividing the rows with the common GCD. */
1017 Vector_Gcd(&MatInv->p[j][0], k, &m1);
1018 Vector_Gcd(&Mat->p[j][0], k, &m2);
1019 value_gcd(gcd, m1, m2);
1020 if (value_notone_p(gcd)) {
1021 for (c = 0; c < k; ++c) {
1022 value_divexact(Mat->p[j][c], Mat->p[j][c], gcd);
1023 value_divexact(MatInv->p[j][c], MatInv->p[j][c], gcd);
1024 }
1025 }
1026 }
1027 }
1028 }
1029
1030 /* Find common denom for each row */
1031 den = (Value *)malloc(k * sizeof(Value));
1032 value_set_si(x, 1);
1033 for (j = 0; j < k; ++j) {
1034 value_init(den[j]);
1035 value_assign(den[j], Mat->p[j][j]);
1036
1037 /* gcd is always positive */
1038 Vector_Gcd(&MatInv->p[j][0], k, &gcd);
1039 value_gcd(gcd, gcd, den[j]);
1040 if (value_neg_p(den[j]))
1041 value_oppose(gcd, gcd); /* make denominator positive */
1042 if (value_notone_p(gcd)) {
1043 for (c = 0; c < k; c++)
1044 value_divexact(MatInv->p[j][c], MatInv->p[j][c], gcd); /* normalize */
1045 value_divexact(den[j], den[j], gcd);
1046 }
1047 value_gcd(gcd, x, den[j]);
1048 value_divexact(m1, den[j], gcd);
1049 value_multiply(x, x, m1);
1050 }
1051 if (value_notone_p(x))
1052 for (j = 0; j < k; ++j) {
1053 for (c = 0; c < k; c++) {
1054 value_division(m1, x, den[j]);
1055 value_multiply(MatInv->p[j][c], MatInv->p[j][c], m1); /* normalize */
1056 }
1057 }
1058
1059 /* Clear all the 'Value' variables */
1060 for (j = 0; j < k; ++j) {
1061 value_clear(den[j]);
1062 }
1063 value_clear(x);
1064 value_clear(gcd);
1065 value_clear(piv);
1066 value_clear(m1);
1067 value_clear(m2);
1068 free(den);
1069
1070 return 1;
1071} /* Matrix_Inverse */
#define value_oppose(ref, val)
Definition: arithmetique.h:552
#define value_notzero_p(val)
Definition: arithmetique.h:576
#define value_divexact(ref, val1, val2)
Definition: arithmetique.h:548
#define value_gcd(ref, val1, val2)
Definition: arithmetique.h:556
#define value_notone_p(val)
Definition: arithmetique.h:578
#define value_one_p(val)
Definition: arithmetique.h:577
#define value_absolute(ref, val)
Definition: arithmetique.h:553
#define value_decrement(ref, val)
Definition: arithmetique.h:546
#define value_zero_p(val)
Definition: arithmetique.h:575
#define value_assign(v1, v2)
Definition: arithmetique.h:482
#define value_set_si(val, i)
Definition: arithmetique.h:483
#define value_addmul(ref, val1, val2)
Definition: arithmetique.h:539
#define value_read(val, str)
Definition: arithmetique.h:486
#define value_clear(val)
Definition: arithmetique.h:485
#define value_division(ref, val1, val2)
Definition: arithmetique.h:547
#define value_multiply(ref, val1, val2)
Definition: arithmetique.h:543
#define value_print(Dst, fmt, val)
Definition: arithmetique.h:487
#define value_lt(v1, v2)
Definition: arithmetique.h:506
#define value_modulus(ref, val1, val2)
Definition: arithmetique.h:549
#define value_subtract(ref, val1, val2)
Definition: arithmetique.h:544
#define value_neg_p(val)
Definition: arithmetique.h:572
#define value_init(val)
Definition: arithmetique.h:481
void errormsg1(const char *f, const char *msgname, const char *msg)
Definition: errormsg.c:25
static int hermite(Matrix *H, Matrix *U, Matrix *Q)
Definition: matrix.c:242
void Vector_Matrix_Product(Value *p1, Matrix *Mat, Value *p2)
Definition: matrix.c:844
void Matrix_Product(Matrix *Mat1, Matrix *Mat2, Matrix *Mat3)
Definition: matrix.c:872
void Matrix_Extend(Matrix *Mat, unsigned NbRows)
Definition: matrix.c:81
void Matrix_Vector_Product(Matrix *Mat, Value *p1, Value *p2)
Definition: matrix.c:812
Matrix * Matrix_Alloc(unsigned NbRows, unsigned NbColumns)
Definition: matrix.c:24
Matrix * Matrix_Read_Input(Matrix *Mat)
Definition: matrix.c:141
int Matrix_Inverse(Matrix *Mat, Matrix *MatInv)
Definition: matrix.c:917
void right_hermite(Matrix *A, Matrix **Hp, Matrix **Up, Matrix **Qp)
Definition: matrix.c:449
Matrix * Matrix_Read(void)
Definition: matrix.c:209
void Matrix_Print(FILE *Dst, const char *Format, Matrix *Mat)
Definition: matrix.c:115
Matrix * Matrix_Read_InputFile(Matrix *Mat, FILE *fp)
Definition: matrix.c:148
void rat_prodmat(Matrix *S, Matrix *X, Matrix *P)
Definition: matrix.c:741
void left_hermite(Matrix *A, Matrix **Hp, Matrix **Qp, Matrix **Up)
Definition: matrix.c:517
int MatInverse(Matrix *Mat, Matrix *MatInv)
Definition: matrix.c:603
Matrix * Matrix_ReadFile(FILE *fp)
Definition: matrix.c:215
void Matrix_Free(Matrix *Mat)
Definition: matrix.c:71
static int nr
Definition: polyparam.c:278
static int n
Definition: polyparam.c:276
char s[128]
Definition: polytest.c:8
Definition: types.h:75
unsigned NbRows
Definition: types.h:76
Value ** p
Definition: types.h:77
int p_Init_size
Definition: types.h:79
unsigned NbColumns
Definition: types.h:76
Value * p_Init
Definition: types.h:78
#define P_VALUE_FMT
Definition: types.h:39
Value * value_alloc(int want, int *got)
Definition: vector.c:781
void Vector_Set(Value *p, int n, unsigned length)
Definition: vector.c:227
void value_free(Value *p, int size)
Definition: vector.c:827
void Vector_Normalize_Positive(Value *p, int length, int pos)
Definition: vector.c:574
void Vector_Gcd(Value *p, unsigned length, Value *min)
Definition: vector.c:485
void Vector_Exchange(Value *p1, Value *p2, unsigned length)
Definition: vector.c:243
void Vector_Copy(Value *p1, Value *p2, unsigned length)
Definition: vector.c:256