Actual source code: dsm.c

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

 15:     /* Local variables */
 16:     PetscInt i,j,maxclq,numgrp;

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

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

 98:     *info = 0;

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

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

110: /*     Determine the degree sequence for the intersection */
111: /*     graph of the columns of A. */

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

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

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

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

128:     if (*maxgrp == *mingrp) {
129:         return(0);
130:     }

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

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

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

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

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

153:         if (*maxgrp == *mingrp) {
154:             return(0);
155:         }
156:     }

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

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

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

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