Actual source code: ido.c

  1: #define PETSCMAT_DLL
  2: /* ido.f -- translated by f2c (version of 25 March 1992  12:58:56).*/

 4:  #include petsc.h
 5:  #include src/mat/color/color.h

  7: static PetscInt c_n1 = -1;

 11: PetscErrorCode MINPACKido(PetscInt *m,PetscInt * n,PetscInt * indrow,PetscInt * jpntr,PetscInt * indcol,PetscInt * ipntr,PetscInt * ndeg,
 12:                PetscInt *list,PetscInt *maxclq, PetscInt *iwa1, PetscInt *iwa2, PetscInt *iwa3, PetscInt *iwa4)
 13: {
 14:     /* System generated locals */
 15:     PetscInt i__1, i__2, i__3, i__4;

 17:     /* Local variables */
 18:     PetscInt jcol = 0, ncomp = 0, ic, ip, jp, ir, maxinc, numinc, numord, maxlst, numwgt, numlst;

 20: /*     Given the sparsity pattern of an m by n matrix A, this */
 21: /*     subroutine determines an incidence-degree ordering of the */
 22: /*     columns of A. */
 23: /*     The incidence-degree ordering is defined for the loopless */
 24: /*     graph G with vertices a(j), j = 1,2,...,n where a(j) is the */
 25: /*     j-th column of A and with edge (a(i),a(j)) if and only if */
 26: /*     columns i and j have a non-zero in the same row position. */
 27: /*     The incidence-degree ordering is determined recursively by */
 28: /*     letting list(k), k = 1,...,n be a column with maximal */
 29: /*     incidence to the subgraph spanned by the ordered columns. */
 30: /*     Among all the columns of maximal incidence, ido chooses a */
 31: /*     column of maximal degree. */
 32: /*     The subroutine statement is */
 33: /*       subroutine ido(m,n,indrow,jpntr,indcol,ipntr,ndeg,list, */
 34: /*                      maxclq,iwa1,iwa2,iwa3,iwa4) */
 35: /*     where */
 36: /*       m is a positive integer input variable set to the number */
 37: /*         of rows of A. */
 38: /*       n is a positive integer input variable set to the number */
 39: /*         of columns of A. */
 40: /*       indrow is an integer input array which contains the row */
 41: /*         indices for the non-zeroes in the matrix A. */
 42: /*       jpntr is an integer input array of length n + 1 which */
 43: /*         specifies the locations of the row indices in indrow. */
 44: /*         The row indices for column j are */
 45: /*               indrow(k), k = jpntr(j),...,jpntr(j+1)-1. */
 46: /*         Note that jpntr(n+1)-1 is then the number of non-zero */
 47: /*         elements of the matrix A. */
 48: /*       indcol is an integer input array which contains the */
 49: /*         column indices for the non-zeroes in the matrix A. */
 50: /*       ipntr is an integer input array of length m + 1 which */
 51: /*         specifies the locations of the column indices in indcol. */
 52: /*         The column indices for row i are */
 53: /*               indcol(k), k = ipntr(i),...,ipntr(i+1)-1. */
 54: /*         Note that ipntr(m+1)-1 is then the number of non-zero */
 55: /*         elements of the matrix A. */
 56: /*       ndeg is an integer input array of length n which specifies */
 57: /*         the degree sequence. The degree of the j-th column */
 58: /*         of A is ndeg(j). */
 59: /*       list is an integer output array of length n which specifies */
 60: /*         the incidence-degree ordering of the columns of A. The j-th */
 61: /*         column in this order is list(j). */
 62: /*       maxclq is an integer output variable set to the size */
 63: /*         of the largest clique found during the ordering. */
 64: /*       iwa1,iwa2,iwa3, and iwa4 are integer work arrays of length n. */
 65: /*     Subprograms called */
 66: /*       MINPACK-supplied ... numsrt */
 67: /*       FORTRAN-supplied ... max */
 68: /*     Argonne National Laboratory. MINPACK Project. August 1984. */
 69: /*     Thomas F. Coleman, Burton S. Garbow, Jorge J. More' */

 71: /*     Sort the degree sequence. */

 74:     /* Parameter adjustments */
 75:     --iwa4;
 76:     --iwa3;
 77:     --iwa2;
 78:     --list;
 79:     --ndeg;
 80:     --ipntr;
 81:     --indcol;
 82:     --jpntr;
 83:     --indrow;

 85:     /* Function Body */
 86:     i__1 = *n - 1;
 87:     MINPACKnumsrt(n, &i__1, &ndeg[1], &c_n1, &iwa4[1], &iwa2[1], &iwa3[1]);

 89: /*     Initialization block. */
 90: /*     Create a doubly-linked list to access the incidences of the */
 91: /*     columns. The pointers for the linked list are as follows. */
 92: /*     Each un-ordered column ic is in a list (the incidence list) */
 93: /*     of columns with the same incidence. */
 94: /*     iwa1(numinc) is the first column in the numinc list */
 95: /*     unless iwa1(numinc) = 0. In this case there are */
 96: /*     no columns in the numinc list. */
 97: /*     iwa2(ic) is the column before ic in the incidence list */
 98: /*     unless iwa2(ic) = 0. In this case ic is the first */
 99: /*     column in this incidence list. */
100: /*     iwa3(ic) is the column after ic in the incidence list */
101: /*     unless iwa3(ic) = 0. In this case ic is the last */
102: /*     column in this incidence list. */
103: /*     If ic is an un-ordered column, then list(ic) is the */
104: /*     incidence of ic to the graph induced by the ordered */
105: /*     columns. If jcol is an ordered column, then list(jcol) */
106: /*     is the incidence-degree order of column jcol. */

108:     maxinc = 0;
109:     for (jp = *n; jp >= 1; --jp) {
110:         ic = iwa4[jp];
111:         iwa1[*n - jp] = 0;
112:         iwa2[ic] = 0;
113:         iwa3[ic] = iwa1[0];
114:         if (iwa1[0] > 0) {
115:             iwa2[iwa1[0]] = ic;
116:         }
117:         iwa1[0] = ic;
118:         iwa4[jp] = 0;
119:         list[jp] = 0;
120:     }

122: /*     Determine the maximal search length for the list */
123: /*     of columns of maximal incidence. */

125:     maxlst = 0;
126:     i__1 = *m;
127:     for (ir = 1; ir <= i__1; ++ir) {
128: /* Computing 2nd power */
129:         i__2 = ipntr[ir + 1] - ipntr[ir];
130:         maxlst += i__2 * i__2;
131:     }
132:     maxlst /= *n;
133:     *maxclq = 0;
134:     numord = 1;

136: /*     Beginning of iteration loop. */

138: L30:

140: /*        Choose a column jcol of maximal degree among the */
141: /*        columns of maximal incidence maxinc. */

143: L40:
144:     jp = iwa1[maxinc];
145:     if (jp > 0) {
146:         goto L50;
147:     }
148:     --maxinc;
149:     goto L40;
150: L50:
151:     numwgt = -1;
152:     i__1 = maxlst;
153:     for (numlst = 1; numlst <= i__1; ++numlst) {
154:         if (ndeg[jp] > numwgt) {
155:             numwgt = ndeg[jp];
156:             jcol = jp;
157:         }
158:         jp = iwa3[jp];
159:         if (jp <= 0) {
160:             goto L70;
161:         }
162:     }
163: L70:
164:     list[jcol] = numord;

166: /*        Update the size of the largest clique */
167: /*        found during the ordering. */

169:     if (!maxinc) {
170:         ncomp = 0;
171:     }
172:     ++ncomp;
173:     if (maxinc + 1 == ncomp) {
174:         *maxclq = PetscMax(*maxclq,ncomp);
175:     }

177: /*        Termination test. */

179:     ++numord;
180:     if (numord > *n) {
181:         goto L100;
182:     }

184: /*        Delete column jcol from the maxinc list. */

186:     if (!iwa2[jcol]) {
187:         iwa1[maxinc] = iwa3[jcol];
188:     } else {
189:         iwa3[iwa2[jcol]] = iwa3[jcol];
190:     }
191:     if (iwa3[jcol] > 0) {
192:         iwa2[iwa3[jcol]] = iwa2[jcol];
193:     }

195: /*        Find all columns adjacent to column jcol. */

197:     iwa4[jcol] = *n;

199: /*        Determine all positions (ir,jcol) which correspond */
200: /*        to non-zeroes in the matrix. */

202:     i__1 = jpntr[jcol + 1] - 1;
203:     for (jp = jpntr[jcol]; jp <= i__1; ++jp) {
204:         ir = indrow[jp];

206: /*           For each row ir, determine all positions (ir,ic) */
207: /*           which correspond to non-zeroes in the matrix. */

209:         i__2 = ipntr[ir + 1] - 1;
210:         for (ip = ipntr[ir]; ip <= i__2; ++ip) {
211:             ic = indcol[ip];

213: /*              Array iwa4 marks columns which are adjacent to */
214: /*              column jcol. */

216:             if (iwa4[ic] < numord) {
217:                 iwa4[ic] = numord;

219: /*                 Update the pointers to the current incidence lists. */

221:                 numinc = list[ic];
222:                 ++list[ic];
223: /* Computing MAX */
224:                 i__3 = maxinc, i__4 = list[ic];
225:                 maxinc = PetscMax(i__3,i__4);

227: /*                 Delete column ic from the numinc list. */

229:                 if (!iwa2[ic]) {
230:                     iwa1[numinc] = iwa3[ic];
231:                 } else {
232:                     iwa3[iwa2[ic]] = iwa3[ic];
233:                 }
234:                 if (iwa3[ic] > 0) {
235:                     iwa2[iwa3[ic]] = iwa2[ic];
236:                 }

238: /*                 Add column ic to the numinc+1 list. */

240:                 iwa2[ic] = 0;
241:                 iwa3[ic] = iwa1[numinc + 1];
242:                 if (iwa1[numinc + 1] > 0) {
243:                     iwa2[iwa1[numinc + 1]] = ic;
244:                 }
245:                 iwa1[numinc + 1] = ic;
246:             }
247:         }
248:     }

250: /*        End of iteration loop. */

252:     goto L30;
253: L100:

255: /*     Invert the array list. */

257:     i__1 = *n;
258:     for (jcol = 1; jcol <= i__1; ++jcol) {
259:         iwa2[list[jcol]] = jcol;
260:     }
261:     i__1 = *n;
262:     for (jp = 1; jp <= i__1; ++jp) {
263:         list[jp] = iwa2[jp];
264:     }
265:     return(0);
266: }