Actual source code: setr.c


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

  4: #include <../src/mat/color/impls/minpack/color.h>

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

 11:   /* Local variables */
 12:   PetscInt jcol, jp, ir;

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

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

 49:   /* Parameter adjustments */
 50:   --iwa;
 51:   --ipntr;
 52:   --indcol;
 53:   --jpntr;
 54:   --indrow;

 56:   /* Function Body */
 57:   i__1 = *m;
 58:   for (ir = 1; ir <= i__1; ++ir) iwa[ir] = 0;

 60:   i__1 = jpntr[*n + 1] - 1;
 61:   for (jp = 1; jp <= i__1; ++jp) ++iwa[indrow[jp]];

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

 65:   ipntr[1] = 1;
 66:   i__1     = *m;

 68:   for (ir = 1; ir <= i__1; ++ir) {
 69:     ipntr[ir + 1] = ipntr[ir] + iwa[ir];
 70:     iwa[ir]       = ipntr[ir];
 71:   }

 73:   /*     Fill indcol. */

 75:   i__1 = *n;
 76:   for (jcol = 1; jcol <= i__1; ++jcol) {
 77:     i__2 = jpntr[jcol + 1] - 1;
 78:     for (jp = jpntr[jcol]; jp <= i__2; ++jp) {
 79:       ir              = indrow[jp];
 80:       indcol[iwa[ir]] = jcol;
 81:       ++iwa[ir];
 82:     }
 83:   }
 84:   return 0;
 85: }