Actual source code: setr.c

  1: #define PETSCMAT_DLL

  3: /* setr.f -- translated by f2c (version of 25 March 1992  12:58:56). */

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

  9: PetscErrorCode MINPACKsetr(PetscInt*m,PetscInt* n,PetscInt* indrow,PetscInt* jpntr,PetscInt* indcol, PetscInt*ipntr,PetscInt* iwa)
 10: {
 11:     /* System generated locals */
 12:     PetscInt i__1, i__2;

 14:     /* Local variables */
 15:     PetscInt jcol, jp, ir;

 17: /*     Given a column-oriented definition of the sparsity pattern */
 18: /*     of an m by n matrix A, this subroutine determines a */
 19: /*     row-oriented definition of the sparsity pattern of A. */
 20: /*     On input the column-oriented definition is specified by */
 21: /*     the arrays indrow and jpntr. On output the row-oriented */
 22: /*     definition is specified by the arrays indcol and ipntr. */
 23: /*     The subroutine statement is */
 24: /*       subroutine setr(m,n,indrow,jpntr,indcol,ipntr,iwa) */
 25: /*     where */
 26: /*       m is a positive integer input variable set to the number */
 27: /*         of rows of A. */
 28: /*       n is a positive integer input variable set to the number */
 29: /*         of columns of A. */
 30: /*       indrow is an integer input array which contains the row */
 31: /*         indices for the non-zeroes in the matrix A. */
 32: /*       jpntr is an integer input array of length n + 1 which */
 33: /*         specifies the locations of the row indices in indrow. */
 34: /*         The row indices for column j are */
 35: /*               indrow(k), k = jpntr(j),...,jpntr(j+1)-1. */
 36: /*         Note that jpntr(n+1)-1 is then the number of non-zero */
 37: /*         elements of the matrix A. */
 38: /*       indcol is an integer output array which contains the */
 39: /*         column indices for the non-zeroes in the matrix A. */
 40: /*       ipntr is an integer output array of length m + 1 which */
 41: /*         specifies the locations of the column indices in indcol. */
 42: /*         The column indices for row i are */
 43: /*               indcol(k), k = ipntr(i),...,ipntr(i+1)-1. */
 44: /*         Note that ipntr(1) is set to 1 and that ipntr(m+1)-1 is */
 45: /*         then the number of non-zero elements of the matrix A. */
 46: /*       iwa is an integer work array of length m. */
 47: /*     Argonne National Laboratory. MINPACK Project. July 1983. */
 48: /*     Thomas F. Coleman, Burton S. Garbow, Jorge J. More' */

 50:     /*     Store in array iwa the counts of non-zeroes in the rows. */

 53:     /* Parameter adjustments */
 54:     --iwa;
 55:     --ipntr;
 56:     --indcol;
 57:     --jpntr;
 58:     --indrow;

 60:     /* Function Body */
 61:     i__1 = *m;
 62:     for (ir = 1; ir <= i__1; ++ir) {
 63:         iwa[ir] = 0;
 64:     }
 65:     i__1 = jpntr[*n + 1] - 1;
 66:     for (jp = 1; jp <= i__1; ++jp) {
 67:         ++iwa[indrow[jp]];
 68:     }

 70:     /*     Set pointers to the start of the rows in indcol. */

 72:     ipntr[1] = 1;
 73:     i__1 = *m;
 74:     for (ir = 1; ir <= i__1; ++ir) {
 75:         ipntr[ir + 1] = ipntr[ir] + iwa[ir];
 76:         iwa[ir] = ipntr[ir];
 77:     }

 79:     /*     Fill indcol. */

 81:     i__1 = *n;
 82:     for (jcol = 1; jcol <= i__1; ++jcol) {
 83:         i__2 = jpntr[jcol + 1] - 1;
 84:         for (jp = jpntr[jcol]; jp <= i__2; ++jp) {
 85:             ir = indrow[jp];
 86:             indcol[iwa[ir]] = jcol;
 87:             ++iwa[ir];
 88:         }
 89:     }
 90:     return(0);
 91: }