Actual source code: dsm.c

  1: /* dsm.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 MINPACKdsm(PetscInt *m, PetscInt *n, PetscInt *npairs, PetscInt *indrow, PetscInt *indcol, PetscInt *ngrp, PetscInt *maxgrp, PetscInt *mingrp, PetscInt *info, PetscInt *ipntr, PetscInt *jpntr, PetscInt *iwa, PetscInt *liwa)
  8: {
  9:   /* System generated locals */
 10:   PetscInt i__1, i__2, i__3;

 12:   /* Local variables */
 13:   PetscInt i, j, maxclq, numgrp;

 15:   /*     Given the sparsity pattern of an m by n matrix A, this */
 16:   /*     subroutine determines a partition of the columns of A */
 17:   /*     consistent with the direct determination of A. */
 18:   /*     The sparsity pattern of the matrix A is specified by */
 19:   /*     the arrays indrow and indcol. On input the indices */
 20:   /*     for the non-zero elements of A are */
 21:   /*           indrow(k),indcol(k), k = 1,2,...,npairs. */
 22:   /*     The (indrow,indcol) pairs may be specified in any order. */
 23:   /*     Duplicate input pairs are permitted, but the subroutine */
 24:   /*     eliminates them. */
 25:   /*     The subroutine partitions the columns of A into groups */
 26:   /*     such that columns in the same group do not have a */
 27:   /*     non-zero in the same row position. A partition of the */
 28:   /*     columns of A with this property is consistent with the */
 29:   /*     direct determination of A. */
 30:   /*     The subroutine statement is */
 31:   /*       subroutine dsm(m,n,npairs,indrow,indcol,ngrp,maxgrp,mingrp, */
 32:   /*                      info,ipntr,jpntr,iwa,liwa) */
 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:   /*       npairs is a positive integer input variable set to the */
 39:   /*         number of (indrow,indcol) pairs used to describe the */
 40:   /*         sparsity pattern of A. */
 41:   /*       indrow is an integer array of length npairs. On input indrow */
 42:   /*         must contain the row indices of the non-zero elements of A. */
 43:   /*         On output indrow is permuted so that the corresponding */
 44:   /*         column indices are in non-decreasing order. The column */
 45:   /*         indices can be recovered from the array jpntr. */
 46:   /*       indcol is an integer array of length npairs. On input indcol */
 47:   /*         must contain the column indices of the non-zero elements of */
 48:   /*         A. On output indcol is permuted so that the corresponding */
 49:   /*         row indices are in non-decreasing order. The row indices */
 50:   /*         can be recovered from the array ipntr. */
 51:   /*       ngrp is an integer output array of length n which specifies */
 52:   /*         the partition of the columns of A. Column jcol belongs */
 53:   /*         to group ngrp(jcol). */
 54:   /*       maxgrp is an integer output variable which specifies the */
 55:   /*         number of groups in the partition of the columns of A. */
 56:   /*       mingrp is an integer output variable which specifies a lower */
 57:   /*         bound for the number of groups in any consistent partition */
 58:   /*         of the columns of A. */
 59:   /*       info is an integer output variable set as follows. For */
 60:   /*         normal termination info = 1. If m, n, or npairs is not */
 61:   /*         positive or liwa is less than max(m,6*n), then info = 0. */
 62:   /*         If the k-th element of indrow is not an integer between */
 63:   /*         1 and m or the k-th element of indcol is not an integer */
 64:   /*         between 1 and n, then info = -k. */
 65:   /*       ipntr is an integer output array of length m + 1 which */
 66:   /*         specifies the locations of the column indices in indcol. */
 67:   /*         The column indices for row i are */
 68:   /*               indcol(k), k = ipntr(i),...,ipntr(i+1)-1. */
 69:   /*         Note that ipntr(m+1)-1 is then the number of non-zero */
 70:   /*         elements of the matrix A. */
 71:   /*       jpntr is an integer output array of length n + 1 which */
 72:   /*         specifies the locations of the row indices in indrow. */
 73:   /*         The row indices for column j are */
 74:   /*               indrow(k), k = jpntr(j),...,jpntr(j+1)-1. */
 75:   /*         Note that jpntr(n+1)-1 is then the number of non-zero */
 76:   /*         elements of the matrix A. */
 77:   /*       iwa is an integer work array of length liwa. */
 78:   /*       liwa is a positive integer input variable not less than */
 79:   /*         max(m,6*n). */
 80:   /*     Subprograms called */
 81:   /*       MINPACK-supplied ... degr,ido,numsrt,seq,setr,slo,srtdat */
 82:   /*       FORTRAN-supplied ... max */
 83:   /*     Argonne National Laboratory. MINPACK Project. December 1984. */
 84:   /*     Thomas F. Coleman, Burton S. Garbow, Jorge J. More' */

 86:   /* Parameter adjustments */
 87:   --iwa;
 88:   --jpntr;
 89:   --ipntr;
 90:   --ngrp;
 91:   --indcol;
 92:   --indrow;

 94:   *info = 0;

 96:   /*     Determine a lower bound for the number of groups. */

 98:   *mingrp = 0;
 99:   i__1    = *m;
100:   for (i = 1; i <= i__1; ++i) {
101:     /* Computing MAX */
102:     i__2    = *mingrp;
103:     i__3    = ipntr[i + 1] - ipntr[i];
104:     *mingrp = PetscMax(i__2, i__3);
105:   }

107:   /*     Determine the degree sequence for the intersection */
108:   /*     graph of the columns of A. */

110:   MINPACKdegr(n, &indrow[1], &jpntr[1], &indcol[1], &ipntr[1], &iwa[*n * 5 + 1], &iwa[*n + 1]);

112:   /*     Color the intersection graph of the columns of A */
113:   /*     with the smallest-last (SL) ordering. */

115:   MINPACKslo(n, &indrow[1], &jpntr[1], &indcol[1], &ipntr[1], &iwa[*n * 5 + 1], &iwa[(*n << 2) + 1], &maxclq, &iwa[1], &iwa[*n + 1], &iwa[(*n << 1) + 1], &iwa[*n * 3 + 1]);
116:   MINPACKseq(n, &indrow[1], &jpntr[1], &indcol[1], &ipntr[1], &iwa[(*n << 2) + 1], &ngrp[1], maxgrp, &iwa[*n + 1]);
117:   *mingrp = PetscMax(*mingrp, maxclq);

119:   /*     Exit if the smallest-last ordering is optimal. */

121:   if (*maxgrp == *mingrp) return 0;

123:   /*     Color the intersection graph of the columns of A */
124:   /*     with the incidence-degree (ID) ordering. */

126:   MINPACKido(m, n, &indrow[1], &jpntr[1], &indcol[1], &ipntr[1], &iwa[*n * 5 + 1], &iwa[(*n << 2) + 1], &maxclq, &iwa[1], &iwa[*n + 1], &iwa[(*n << 1) + 1], &iwa[*n * 3 + 1]);
127:   MINPACKseq(n, &indrow[1], &jpntr[1], &indcol[1], &ipntr[1], &iwa[(*n << 2) + 1], &iwa[1], &numgrp, &iwa[*n + 1]);
128:   *mingrp = PetscMax(*mingrp, maxclq);

130:   /*     Retain the better of the two orderings so far. */

132:   if (numgrp < *maxgrp) {
133:     *maxgrp = numgrp;
134:     i__1    = *n;
135:     for (j = 1; j <= i__1; ++j) ngrp[j] = iwa[j];

137:     /*        Exit if the incidence-degree ordering is optimal. */

139:     if (*maxgrp == *mingrp) return 0;
140:   }

142:   /*     Color the intersection graph of the columns of A */
143:   /*     with the largest-first (LF) ordering. */

145:   i__1 = *n - 1;
146:   MINPACKnumsrt(n, &i__1, &iwa[*n * 5 + 1], &c_n1, &iwa[(*n << 2) + 1], &iwa[(*n << 1) + 1], &iwa[*n + 1]);
147:   MINPACKseq(n, &indrow[1], &jpntr[1], &indcol[1], &ipntr[1], &iwa[(*n << 2) + 1], &iwa[1], &numgrp, &iwa[*n + 1]);

149:   /*     Retain the best of the three orderings and exit. */

151:   if (numgrp < *maxgrp) {
152:     *maxgrp = numgrp;
153:     i__1    = *n;
154:     for (j = 1; j <= i__1; ++j) ngrp[j] = iwa[j];
155:   }
156:   return 0;
157: }