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