Actual source code: ido.c
1: /* ido.f -- translated by f2c (version of 25 March 1992 12:58:56).*/
3: #include <../src/mat/color/impls/minpack/color.h>
5: static PetscInt c_n1 = -1;
7: PetscErrorCode MINPACKido(PetscInt *m, PetscInt *n, const PetscInt *indrow, const PetscInt *jpntr, const PetscInt *indcol, const PetscInt *ipntr, PetscInt *ndeg, PetscInt *list, PetscInt *maxclq, PetscInt *iwa1, PetscInt *iwa2, PetscInt *iwa3, PetscInt *iwa4)
8: {
9: /* System generated locals */
10: PetscInt i__1, i__2, i__3, i__4;
12: /* Local variables */
13: PetscInt jcol = 0, ncomp = 0, ic, ip, jp, ir, maxinc, numinc, numord, maxlst, numwgt, numlst;
15: /* Given the sparsity pattern of an m by n matrix A, this */
16: /* subroutine determines an incidence-degree ordering of the */
17: /* columns of A. */
18: /* The incidence-degree ordering is defined for the loopless */
19: /* graph G with vertices a(j), j = 1,2,...,n where a(j) is the */
20: /* j-th column of A and with edge (a(i),a(j)) if and only if */
21: /* columns i and j have a non-zero in the same row position. */
22: /* The incidence-degree ordering is determined recursively by */
23: /* letting list(k), k = 1,...,n be a column with maximal */
24: /* incidence to the subgraph spanned by the ordered columns. */
25: /* Among all the columns of maximal incidence, ido chooses a */
26: /* column of maximal degree. */
27: /* The subroutine statement is */
28: /* subroutine ido(m,n,indrow,jpntr,indcol,ipntr,ndeg,list, */
29: /* maxclq,iwa1,iwa2,iwa3,iwa4) */
30: /* where */
31: /* m is a positive integer input variable set to the number */
32: /* of rows of A. */
33: /* n is a positive integer input variable set to the number */
34: /* of columns of A. */
35: /* indrow is an integer input array which contains the row */
36: /* indices for the non-zeroes in the matrix A. */
37: /* jpntr is an integer input array of length n + 1 which */
38: /* specifies the locations of the row indices in indrow. */
39: /* The row indices for column j are */
40: /* indrow(k), k = jpntr(j),...,jpntr(j+1)-1. */
41: /* Note that jpntr(n+1)-1 is then the number of non-zero */
42: /* elements of the matrix A. */
43: /* indcol is an integer input array which contains the */
44: /* column indices for the non-zeroes in the matrix A. */
45: /* ipntr is an integer input array of length m + 1 which */
46: /* specifies the locations of the column indices in indcol. */
47: /* The column indices for row i are */
48: /* indcol(k), k = ipntr(i),...,ipntr(i+1)-1. */
49: /* Note that ipntr(m+1)-1 is then the number of non-zero */
50: /* elements of the matrix A. */
51: /* ndeg is an integer input array of length n which specifies */
52: /* the degree sequence. The degree of the j-th column */
53: /* of A is ndeg(j). */
54: /* list is an integer output array of length n which specifies */
55: /* the incidence-degree ordering of the columns of A. The j-th */
56: /* column in this order is list(j). */
57: /* maxclq is an integer output variable set to the size */
58: /* of the largest clique found during the ordering. */
59: /* iwa1,iwa2,iwa3, and iwa4 are integer work arrays of length n. */
60: /* Subprograms called */
61: /* MINPACK-supplied ... numsrt */
62: /* FORTRAN-supplied ... max */
63: /* Argonne National Laboratory. MINPACK Project. August 1984. */
64: /* Thomas F. Coleman, Burton S. Garbow, Jorge J. More' */
66: /* Sort the degree sequence. */
68: /* Parameter adjustments */
69: --iwa4;
70: --iwa3;
71: --iwa2;
72: --list;
73: --ndeg;
74: --ipntr;
75: --indcol;
76: --jpntr;
77: --indrow;
79: /* Function Body */
80: i__1 = *n - 1;
81: MINPACKnumsrt(n, &i__1, &ndeg[1], &c_n1, &iwa4[1], &iwa2[1], &iwa3[1]);
83: /* Initialization block. */
84: /* Create a doubly-linked list to access the incidences of the */
85: /* columns. The pointers for the linked list are as follows. */
86: /* Each un-ordered column ic is in a list (the incidence list) */
87: /* of columns with the same incidence. */
88: /* iwa1(numinc) is the first column in the numinc list */
89: /* unless iwa1(numinc) = 0. In this case there are */
90: /* no columns in the numinc list. */
91: /* iwa2(ic) is the column before ic in the incidence list */
92: /* unless iwa2(ic) = 0. In this case ic is the first */
93: /* column in this incidence list. */
94: /* iwa3(ic) is the column after ic in the incidence list */
95: /* unless iwa3(ic) = 0. In this case ic is the last */
96: /* column in this incidence list. */
97: /* If ic is an un-ordered column, then list(ic) is the */
98: /* incidence of ic to the graph induced by the ordered */
99: /* columns. If jcol is an ordered column, then list(jcol) */
100: /* is the incidence-degree order of column jcol. */
102: maxinc = 0;
103: for (jp = *n; jp >= 1; --jp) {
104: ic = iwa4[jp];
105: iwa1[*n - jp] = 0;
106: iwa2[ic] = 0;
107: iwa3[ic] = iwa1[0];
108: if (iwa1[0] > 0) iwa2[iwa1[0]] = ic;
109: iwa1[0] = ic;
110: iwa4[jp] = 0;
111: list[jp] = 0;
112: }
114: /* Determine the maximal search length for the list */
115: /* of columns of maximal incidence. */
117: maxlst = 0;
118: i__1 = *m;
119: for (ir = 1; ir <= i__1; ++ir) {
120: /* Computing 2nd power */
121: i__2 = ipntr[ir + 1] - ipntr[ir];
122: maxlst += i__2 * i__2;
123: }
124: maxlst /= *n;
125: *maxclq = 0;
126: numord = 1;
128: /* Beginning of iteration loop. */
130: L30:
132: /* Choose a column jcol of maximal degree among the */
133: /* columns of maximal incidence maxinc. */
135: L40:
136: jp = iwa1[maxinc];
137: if (jp > 0) goto L50;
138: --maxinc;
139: goto L40;
140: L50:
141: numwgt = -1;
142: i__1 = maxlst;
143: for (numlst = 1; numlst <= i__1; ++numlst) {
144: if (ndeg[jp] > numwgt) {
145: numwgt = ndeg[jp];
146: jcol = jp;
147: }
148: jp = iwa3[jp];
149: if (jp <= 0) goto L70;
150: }
151: L70:
152: list[jcol] = numord;
154: /* Update the size of the largest clique */
155: /* found during the ordering. */
157: if (!maxinc) ncomp = 0;
158: ++ncomp;
159: if (maxinc + 1 == ncomp) *maxclq = PetscMax(*maxclq, ncomp);
161: /* Termination test. */
163: ++numord;
164: if (numord > *n) goto L100;
166: /* Delete column jcol from the maxinc list. */
168: if (!iwa2[jcol]) iwa1[maxinc] = iwa3[jcol];
169: else iwa3[iwa2[jcol]] = iwa3[jcol];
171: if (iwa3[jcol] > 0) iwa2[iwa3[jcol]] = iwa2[jcol];
173: /* Find all columns adjacent to column jcol. */
175: iwa4[jcol] = *n;
177: /* Determine all positions (ir,jcol) which correspond */
178: /* to non-zeroes in the matrix. */
180: i__1 = jpntr[jcol + 1] - 1;
181: for (jp = jpntr[jcol]; jp <= i__1; ++jp) {
182: ir = indrow[jp];
184: /* For each row ir, determine all positions (ir,ic) */
185: /* which correspond to non-zeroes in the matrix. */
187: i__2 = ipntr[ir + 1] - 1;
188: for (ip = ipntr[ir]; ip <= i__2; ++ip) {
189: ic = indcol[ip];
191: /* Array iwa4 marks columns which are adjacent to */
192: /* column jcol. */
194: if (iwa4[ic] < numord) {
195: iwa4[ic] = numord;
197: /* Update the pointers to the current incidence lists. */
199: numinc = list[ic];
200: ++list[ic];
201: /* Computing MAX */
202: i__3 = maxinc;
203: i__4 = list[ic];
204: maxinc = PetscMax(i__3, i__4);
206: /* Delete column ic from the numinc list. */
208: if (!iwa2[ic]) iwa1[numinc] = iwa3[ic];
209: else iwa3[iwa2[ic]] = iwa3[ic];
211: if (iwa3[ic] > 0) iwa2[iwa3[ic]] = iwa2[ic];
213: /* Add column ic to the numinc+1 list. */
215: iwa2[ic] = 0;
216: iwa3[ic] = iwa1[numinc + 1];
217: if (iwa1[numinc + 1] > 0) iwa2[iwa1[numinc + 1]] = ic;
218: iwa1[numinc + 1] = ic;
219: }
220: }
221: }
223: /* End of iteration loop. */
225: goto L30;
226: L100:
228: /* Invert the array list. */
230: i__1 = *n;
231: for (jcol = 1; jcol <= i__1; ++jcol) iwa2[list[jcol]] = jcol;
232: i__1 = *n;
233: for (jp = 1; jp <= i__1; ++jp) list[jp] = iwa2[jp];
234: return 0;
235: }