Actual source code: dsm.c

  1: #define PETSCMAT_DLL
  2: /* dsm.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 MINPACKdsm(PetscInt *m,PetscInt *n,PetscInt *npairs,PetscInt *indrow,PetscInt *indcol,PetscInt *ngrp,PetscInt *maxgrp,
 12:                           PetscInt *mingrp,PetscInt *info,PetscInt *ipntr,PetscInt *jpntr,PetscInt *iwa,PetscInt *liwa)
 13: {
 14:     /* System generated locals */
 15:     PetscInt i__1,i__2,i__3;

 17:     /* Local variables */
 18:     PetscInt i,j,maxclq,numgrp;

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

 92:     /* Parameter adjustments */
 93:     --iwa;
 94:     --jpntr;
 95:     --ipntr;
 96:     --ngrp;
 97:     --indcol;
 98:     --indrow;

100:     *info = 0;

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

104:     *mingrp = 0;
105:     i__1 = *m;
106:     for (i = 1; i <= i__1; ++i) {
107: /* Computing MAX */
108:         i__2 = *mingrp,i__3 = ipntr[i + 1] - ipntr[i];
109:         *mingrp = PetscMax(i__2,i__3);
110:     }

112: /*     Determine the degree sequence for the intersection */
113: /*     graph of the columns of A. */

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

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

121:     MINPACKslo(n,&indrow[1],&jpntr[1],&indcol[1],&ipntr[1],&iwa[*n * 5 + 1],&
122:             iwa[(*n << 2) + 1],&maxclq,&iwa[1],&iwa[*n + 1],&iwa[(*n << 1)
123:              + 1],&iwa[*n * 3 + 1]);
124:     MINPACKseq(n,&indrow[1],&jpntr[1],&indcol[1],&ipntr[1],&iwa[(*n << 2) + 1],
125:              &ngrp[1],maxgrp,&iwa[*n + 1]);
126:     *mingrp = PetscMax(*mingrp,maxclq);

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

130:     if (*maxgrp == *mingrp) {
131:         return(0);
132:     }

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

137:     MINPACKido(m,n,&indrow[1],&jpntr[1],&indcol[1],&ipntr[1],&iwa[*n * 5 + 1],
138:              &iwa[(*n << 2) + 1],&maxclq,&iwa[1],&iwa[*n + 1],&iwa[(*n <<
139:             1) + 1],&iwa[*n * 3 + 1]);
140:     MINPACKseq(n,&indrow[1],&jpntr[1],&indcol[1],&ipntr[1],&iwa[(*n << 2) + 1],
141:              &iwa[1],&numgrp,&iwa[*n + 1]);
142:     *mingrp = PetscMax(*mingrp,maxclq);

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

146:     if (numgrp < *maxgrp) {
147:         *maxgrp = numgrp;
148:         i__1 = *n;
149:         for (j = 1; j <= i__1; ++j) {
150:             ngrp[j] = iwa[j];
151:         }

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

155:         if (*maxgrp == *mingrp) {
156:             return(0);
157:         }
158:     }

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

163:     i__1 = *n - 1;
164:     MINPACKnumsrt(n,&i__1,&iwa[*n * 5 + 1],&c_n1,&iwa[(*n << 2) + 1],&iwa[(*n
165:             << 1) + 1],&iwa[*n + 1]);
166:     MINPACKseq(n,&indrow[1],&jpntr[1],&indcol[1],&ipntr[1],&iwa[(*n << 2) + 1],
167:              &iwa[1],&numgrp,&iwa[*n + 1]);

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

171:     if (numgrp < *maxgrp) {
172:         *maxgrp = numgrp;
173:         i__1 = *n;
174:         for (j = 1; j <= i__1; ++j) {
175:             ngrp[j] = iwa[j];
176:         }
177:     }
178:     return(0);
179: }