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