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