Actual source code: degr.c
2: /* degr.f -- translated by f2c (version of 25 March 1992 12:58:56). */
4: #include <../src/mat/color/impls/minpack/color.h>
6: PetscErrorCode MINPACKdegr(PetscInt *n, const PetscInt *indrow, const PetscInt *jpntr, const PetscInt *indcol, const PetscInt *ipntr, PetscInt *ndeg, PetscInt *iwa)
7: {
8: /* System generated locals */
9: PetscInt i__1, i__2, i__3;
11: /* Local variables */
12: PetscInt jcol, ic, ip, jp, ir;
14: /* subroutine degr */
15: /* Given the sparsity pattern of an m by n matrix A, */
16: /* this subroutine determines the degree sequence for */
17: /* the intersection graph of the columns of A. */
18: /* In graph-theory terminology, the intersection graph of */
19: /* the columns of A is the loopless graph G with vertices */
20: /* a(j), j = 1,2,...,n where a(j) is the j-th column of A */
21: /* and with edge (a(i),a(j)) if and only if columns i and j */
22: /* have a non-zero in the same row position. */
23: /* Note that the value of m is not needed by degr and is */
24: /* therefore not present in the subroutine statement. */
25: /* The subroutine statement is */
26: /* subroutine degr(n,indrow,jpntr,indcol,ipntr,ndeg,iwa) */
27: /* where */
28: /* n is a positive integer input variable set to the number */
29: /* of columns of A. */
30: /* indrow is an integer input array which contains the row */
31: /* indices for the non-zeroes in the matrix A. */
32: /* jpntr is an integer input array of length n + 1 which */
33: /* specifies the locations of the row indices in indrow. */
34: /* The row indices for column j are */
35: /* indrow(k), k = jpntr(j),...,jpntr(j+1)-1. */
36: /* Note that jpntr(n+1)-1 is then the number of non-zero */
37: /* elements of the matrix A. */
38: /* indcol is an integer input array which contains the */
39: /* column indices for the non-zeroes in the matrix A. */
40: /* ipntr is an integer input array of length m + 1 which */
41: /* specifies the locations of the column indices in indcol. */
42: /* The column indices for row i are */
43: /* indcol(k), k = ipntr(i),...,ipntr(i+1)-1. */
44: /* Note that ipntr(m+1)-1 is then the number of non-zero */
45: /* elements of the matrix A. */
46: /* ndeg is an integer output array of length n which */
47: /* specifies the degree sequence. The degree of the */
48: /* j-th column of A is ndeg(j). */
49: /* iwa is an integer work array of length n. */
50: /* Argonne National Laboratory. MINPACK Project. July 1983. */
51: /* Thomas F. Coleman, Burton S. Garbow, Jorge J. More' */
53: /* Parameter adjustments */
54: --iwa;
55: --ndeg;
56: --ipntr;
57: --indcol;
58: --jpntr;
59: --indrow;
61: /* Function Body */
62: i__1 = *n;
63: for (jp = 1; jp <= i__1; ++jp) {
64: ndeg[jp] = 0;
65: iwa[jp] = 0;
66: }
68: /* Compute the degree sequence by determining the contributions */
69: /* to the degrees from the current(jcol) column and further */
70: /* columns which have not yet been considered. */
72: i__1 = *n;
73: for (jcol = 2; jcol <= i__1; ++jcol) {
74: iwa[jcol] = *n;
76: /* Determine all positions (ir,jcol) which correspond */
77: /* to non-zeroes in the matrix. */
79: i__2 = jpntr[jcol + 1] - 1;
80: for (jp = jpntr[jcol]; jp <= i__2; ++jp) {
81: ir = indrow[jp];
83: /* For each row ir, determine all positions (ir,ic) */
84: /* which correspond to non-zeroes in the matrix. */
86: i__3 = ipntr[ir + 1] - 1;
87: for (ip = ipntr[ir]; ip <= i__3; ++ip) {
88: ic = indcol[ip];
90: /* Array iwa marks columns which have contributed to */
91: /* the degree count of column jcol. Update the degree */
92: /* counts of these columns as well as column jcol. */
94: if (iwa[ic] < jcol) {
95: iwa[ic] = jcol;
96: ++ndeg[ic];
97: ++ndeg[jcol];
98: }
99: }
100: }
101: }
102: return 0;
103: }